xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision eec179cf895b6fbcd6a0b58694319392b06e361b)
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 
14d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C, PetscInt imax, IS is[], PetscInt ov)
15d71ae5a4SJacob Faibussowitsch {
16d0f46423SBarry Smith   PetscInt i, N = C->cmap->N, bs = C->rmap->bs;
1736f4e84dSSatish Balay   IS      *is_new;
1836f4e84dSSatish Balay 
193a40ed3dSBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(imax, &is_new));
21df36731bSBarry Smith   /* Convert the indices into block format */
229566063dSJacob Faibussowitsch   PetscCall(ISCompressIndicesGeneral(N, C->rmap->n, bs, imax, is, is_new));
2308401ef6SPierre Jolivet   PetscCheck(ov >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified");
2448a46eb9SPierre Jolivet   for (i = 0; i < ov; ++i) PetscCall(MatIncreaseOverlap_MPIBAIJ_Once(C, imax, is_new));
259566063dSJacob Faibussowitsch   for (i = 0; i < imax; i++) PetscCall(ISDestroy(&is[i]));
269566063dSJacob Faibussowitsch   PetscCall(ISExpandIndicesGeneral(N, N, bs, imax, is_new, is));
279566063dSJacob Faibussowitsch   for (i = 0; i < imax; i++) PetscCall(ISDestroy(&is_new[i]));
289566063dSJacob Faibussowitsch   PetscCall(PetscFree(is_new));
293a40ed3dSBarry Smith   PetscFunctionReturn(0);
30d5b485abSSatish Balay }
31d5b485abSSatish Balay 
32d5b485abSSatish Balay /*
33d5b485abSSatish Balay   Sample message format:
34d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
350e9b0e7eSHong Zhang   to index sets is[1], is[5]
36d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
37d5b485abSSatish Balay   -----------
38d5b485abSSatish Balay   mesg [1] = 1 => is[1]
39d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
40d5b485abSSatish Balay   -----------
41d5b485abSSatish Balay   mesg [5] = 5  => is[5]
42d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
43d5b485abSSatish Balay   -----------
44d5b485abSSatish Balay   mesg [7]
450e9b0e7eSHong Zhang   mesg [n]  data(is[1])
46d5b485abSSatish Balay   -----------
47d5b485abSSatish Balay   mesg[n+1]
48d5b485abSSatish Balay   mesg[m]  data(is[5])
49d5b485abSSatish Balay   -----------
50d5b485abSSatish Balay 
51d5b485abSSatish Balay   Notes:
52d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
532d4ee042Sprj-   nrqr - no of requests received (which have to be or which have been processed)
54d5b485abSSatish Balay */
55d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C, PetscInt imax, IS is[])
56d71ae5a4SJacob Faibussowitsch {
57df36731bSBarry Smith   Mat_MPIBAIJ     *c = (Mat_MPIBAIJ *)C->data;
585d0c19d7SBarry Smith   const PetscInt **idx, *idx_i;
5924c049a4SHong Zhang   PetscInt        *n, *w3, *w4, **data, len;
60b24ad042SBarry Smith   PetscMPIInt      size, rank, tag1, tag2, *w2, *w1, nrqr;
61131c27b5Sprj-   PetscInt         Mbs, i, j, k, **rbuf, row, nrqs, msz, **outdat, **ptr;
628f9f447aSHong Zhang   PetscInt        *ctr, *pa, *tmp, *isz, *isz1, **xdata, **rbuf2, *d_p;
63131c27b5Sprj-   PetscMPIInt     *onodes1, *olengths1, *onodes2, *olengths2, proc = -1;
64f1af5d2fSBarry Smith   PetscBT         *table;
65749130a8SStefano Zampini   MPI_Comm         comm, *iscomms;
66d5b485abSSatish Balay   MPI_Request     *s_waits1, *r_waits1, *s_waits2, *r_waits2;
678f9f447aSHong Zhang   char            *t_p;
68d5b485abSSatish Balay 
693a40ed3dSBarry Smith   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
71d5b485abSSatish Balay   size = c->size;
72d5b485abSSatish Balay   rank = c->rank;
73df36731bSBarry Smith   Mbs  = c->Mbs;
74d5b485abSSatish Balay 
759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag1));
769566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2));
77c7dd2924SSatish Balay 
789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(imax + 1, (PetscInt ***)&idx, imax, &n));
7924c049a4SHong Zhang 
80d5b485abSSatish Balay   for (i = 0; i < imax; i++) {
819566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[i], &idx[i]));
829566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(is[i], &n[i]));
83d5b485abSSatish Balay   }
84d5b485abSSatish Balay 
85d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
86d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
879566063dSJacob Faibussowitsch   PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4));
88d5b485abSSatish Balay   for (i = 0; i < imax; i++) {
899566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(w4, size)); /* initialise work vector*/
90d5b485abSSatish Balay     idx_i = idx[i];
91d5b485abSSatish Balay     len   = n[i];
92d5b485abSSatish Balay     for (j = 0; j < len; j++) {
93d5b485abSSatish Balay       row = idx_i[j];
9408401ef6SPierre Jolivet       PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index set cannot have negative entries");
959566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc));
96d5b485abSSatish Balay       w4[proc]++;
97d5b485abSSatish Balay     }
98d5b485abSSatish Balay     for (j = 0; j < size; j++) {
999371c9d4SSatish Balay       if (w4[j]) {
1009371c9d4SSatish Balay         w1[j] += w4[j];
1019371c9d4SSatish Balay         w3[j]++;
1029371c9d4SSatish Balay       }
103d5b485abSSatish Balay     }
104d5b485abSSatish Balay   }
105d5b485abSSatish Balay 
106d5b485abSSatish Balay   nrqs     = 0; /* no of outgoing messages */
107d5b485abSSatish Balay   msz      = 0; /* total mesg length (for all proc */
1080e9b0e7eSHong Zhang   w1[rank] = 0; /* no mesg sent to itself */
109d5b485abSSatish Balay   w3[rank] = 0;
110d5b485abSSatish Balay   for (i = 0; i < size; i++) {
1119371c9d4SSatish Balay     if (w1[i]) {
1129371c9d4SSatish Balay       w2[i] = 1;
1139371c9d4SSatish Balay       nrqs++;
1149371c9d4SSatish Balay     } /* there exists a message to proc i */
115d5b485abSSatish Balay   }
116d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
1179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqs, &pa));
118d5b485abSSatish Balay   for (i = 0, j = 0; i < size; i++) {
1199371c9d4SSatish Balay     if (w1[i]) {
1209371c9d4SSatish Balay       pa[j] = i;
1219371c9d4SSatish Balay       j++;
1229371c9d4SSatish Balay     }
123d5b485abSSatish Balay   }
124d5b485abSSatish Balay 
125d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
126d5b485abSSatish Balay   for (i = 0; i < nrqs; i++) {
127d5b485abSSatish Balay     j = pa[i];
128d5b485abSSatish Balay     w1[j] += w2[j] + 2 * w3[j];
129d5b485abSSatish Balay     msz += w1[j];
130d5b485abSSatish Balay   }
131d5b485abSSatish Balay 
132c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
1339566063dSJacob Faibussowitsch   PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr));
1349566063dSJacob Faibussowitsch   PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1));
135d5b485abSSatish Balay 
136c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
1379566063dSJacob Faibussowitsch   PetscCall(PetscPostIrecvInt(comm, tag1, nrqr, onodes1, olengths1, &rbuf, &r_waits1));
138d5b485abSSatish Balay 
139d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
1409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(size, &outdat, size, &ptr, msz, &tmp, size, &ctr));
1419566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(outdat, size));
1429566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(ptr, size));
143d5b485abSSatish Balay   {
144b24ad042SBarry Smith     PetscInt *iptr = tmp, ict = 0;
145d5b485abSSatish Balay     for (i = 0; i < nrqs; i++) {
146d5b485abSSatish Balay       j = pa[i];
147d5b485abSSatish Balay       iptr += ict;
148d5b485abSSatish Balay       outdat[j] = iptr;
149d5b485abSSatish Balay       ict       = w1[j];
150d5b485abSSatish Balay     }
151d5b485abSSatish Balay   }
152d5b485abSSatish Balay 
153d5b485abSSatish Balay   /* Form the outgoing messages */
154d5b485abSSatish Balay   /*plug in the headers*/
155d5b485abSSatish Balay   for (i = 0; i < nrqs; i++) {
156d5b485abSSatish Balay     j            = pa[i];
157d5b485abSSatish Balay     outdat[j][0] = 0;
1589566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(outdat[j] + 1, 2 * w3[j]));
159d5b485abSSatish Balay     ptr[j] = outdat[j] + 2 * w3[j] + 1;
160d5b485abSSatish Balay   }
161d5b485abSSatish Balay 
162d5b485abSSatish Balay   /* Memory for doing local proc's work*/
163d5b485abSSatish Balay   {
1649566063dSJacob Faibussowitsch     PetscCall(PetscCalloc5(imax, &table, imax, &data, imax, &isz, Mbs * imax, &d_p, (Mbs / PETSC_BITS_PER_BYTE + 1) * imax, &t_p));
165d5b485abSSatish Balay 
166bbd702dbSSatish Balay     for (i = 0; i < imax; i++) {
167b6410449SSatish Balay       table[i] = t_p + (Mbs / PETSC_BITS_PER_BYTE + 1) * i;
168bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
169d5b485abSSatish Balay     }
170d5b485abSSatish Balay   }
171d5b485abSSatish Balay 
172d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
173d5b485abSSatish Balay   {
174b24ad042SBarry Smith     PetscInt n_i, *data_i, isz_i, *outdat_j, ctr_j;
175f1af5d2fSBarry Smith     PetscBT  table_i;
176d5b485abSSatish Balay 
177d5b485abSSatish Balay     for (i = 0; i < imax; i++) {
1789566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctr, size));
179d5b485abSSatish Balay       n_i     = n[i];
180d5b485abSSatish Balay       table_i = table[i];
181d5b485abSSatish Balay       idx_i   = idx[i];
182d5b485abSSatish Balay       data_i  = data[i];
183d5b485abSSatish Balay       isz_i   = isz[i];
184d5b485abSSatish Balay       for (j = 0; j < n_i; j++) { /* parse the indices of each IS */
185d5b485abSSatish Balay         row = idx_i[j];
1869566063dSJacob Faibussowitsch         PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc));
187d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
188d5b485abSSatish Balay           ctr[proc]++;
189d5b485abSSatish Balay           *ptr[proc] = row;
190d5b485abSSatish Balay           ptr[proc]++;
191d6b45a43SBarry Smith         } else { /* Update the local table */
19226fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row;
193d5b485abSSatish Balay         }
194d5b485abSSatish Balay       }
195d5b485abSSatish Balay       /* Update the headers for the current IS */
196d5b485abSSatish Balay       for (j = 0; j < size; j++) { /* Can Optimise this loop by using pa[] */
197d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
198d5b485abSSatish Balay           outdat_j            = outdat[j];
199d5b485abSSatish Balay           k                   = ++outdat_j[0];
200d5b485abSSatish Balay           outdat_j[2 * k]     = ctr_j;
201d5b485abSSatish Balay           outdat_j[2 * k - 1] = i;
202d5b485abSSatish Balay         }
203d5b485abSSatish Balay       }
204d5b485abSSatish Balay       isz[i] = isz_i;
205d5b485abSSatish Balay     }
206d5b485abSSatish Balay   }
207d5b485abSSatish Balay 
208d5b485abSSatish Balay   /*  Now  post the sends */
2099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqs, &s_waits1));
210d5b485abSSatish Balay   for (i = 0; i < nrqs; ++i) {
211d5b485abSSatish Balay     j = pa[i];
2129566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Isend(outdat[j], w1[j], MPIU_INT, j, tag1, comm, s_waits1 + i));
213d5b485abSSatish Balay   }
214d5b485abSSatish Balay 
215d5b485abSSatish Balay   /* No longer need the original indices*/
21648a46eb9SPierre Jolivet   for (i = 0; i < imax; ++i) PetscCall(ISRestoreIndices(is[i], idx + i));
2179566063dSJacob Faibussowitsch   PetscCall(PetscFree2(*(PetscInt ***)&idx, n));
218d5b485abSSatish Balay 
2199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(imax, &iscomms));
220d5b485abSSatish Balay   for (i = 0; i < imax; ++i) {
2219566063dSJacob Faibussowitsch     PetscCall(PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]), &iscomms[i], NULL));
2229566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is[i]));
223d5b485abSSatish Balay   }
224d5b485abSSatish Balay 
225d5b485abSSatish Balay   /* Do Local work*/
2269566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap_MPIBAIJ_Local(C, imax, table, isz, data));
227d5b485abSSatish Balay 
228d5b485abSSatish Balay   /* Receive messages*/
2299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE));
2309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE));
231d5b485abSSatish Balay 
232d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
2339566063dSJacob Faibussowitsch   PetscCall(PetscFree4(outdat, ptr, tmp, ctr));
2349566063dSJacob Faibussowitsch   PetscCall(PetscFree4(w1, w2, w3, w4));
235d5b485abSSatish Balay 
2369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqr, &xdata));
2379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqr, &isz1));
2389566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap_MPIBAIJ_Receive(C, nrqr, rbuf, xdata, isz1));
239175c2bc5SJunchao Zhang   if (rbuf) {
2409566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf[0]));
2419566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf));
242175c2bc5SJunchao Zhang   }
243d5b485abSSatish Balay 
244d5b485abSSatish Balay   /* Send the data back*/
245d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
246d5b485abSSatish Balay   {
247b24ad042SBarry Smith     PetscMPIInt *rw1;
248d5b485abSSatish Balay 
2499566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(size, &rw1));
250c7dd2924SSatish Balay 
251d5b485abSSatish Balay     for (i = 0; i < nrqr; ++i) {
252175c2bc5SJunchao Zhang       proc      = onodes1[i];
253d5b485abSSatish Balay       rw1[proc] = isz1[i];
254d5b485abSSatish Balay     }
255d5b485abSSatish Balay 
256c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
2579566063dSJacob Faibussowitsch     PetscCall(PetscGatherMessageLengths(comm, nrqr, nrqs, rw1, &onodes2, &olengths2));
2589566063dSJacob Faibussowitsch     PetscCall(PetscFree(rw1));
259c7dd2924SSatish Balay   }
260c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
2619566063dSJacob Faibussowitsch   PetscCall(PetscPostIrecvInt(comm, tag2, nrqs, onodes2, olengths2, &rbuf2, &r_waits2));
262d5b485abSSatish Balay 
263d5b485abSSatish Balay   /*  Now  post the sends */
2649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqr, &s_waits2));
265d5b485abSSatish Balay   for (i = 0; i < nrqr; ++i) {
266175c2bc5SJunchao Zhang     j = onodes1[i];
2679566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Isend(xdata[i], isz1[i], MPIU_INT, j, tag2, comm, s_waits2 + i));
268d5b485abSSatish Balay   }
269d5b485abSSatish Balay 
2709566063dSJacob Faibussowitsch   PetscCall(PetscFree(onodes1));
2719566063dSJacob Faibussowitsch   PetscCall(PetscFree(olengths1));
272175c2bc5SJunchao Zhang 
273d5b485abSSatish Balay   /* receive work done on other processors*/
274d5b485abSSatish Balay   {
275b24ad042SBarry Smith     PetscMPIInt idex;
276b24ad042SBarry Smith     PetscInt    is_no, ct1, max, *rbuf2_i, isz_i, *data_i, jmax;
277f1af5d2fSBarry Smith     PetscBT     table_i;
278d5b485abSSatish Balay 
279d5b485abSSatish Balay     for (i = 0; i < nrqs; ++i) {
2809566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Waitany(nrqs, r_waits2, &idex, MPI_STATUS_IGNORE));
281d5b485abSSatish Balay       /* Process the message*/
282999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
283d5b485abSSatish Balay       ct1     = 2 * rbuf2_i[0] + 1;
284999d9058SBarry Smith       jmax    = rbuf2[idex][0];
285d5b485abSSatish Balay       for (j = 1; j <= jmax; j++) {
286d5b485abSSatish Balay         max     = rbuf2_i[2 * j];
287d5b485abSSatish Balay         is_no   = rbuf2_i[2 * j - 1];
288d5b485abSSatish Balay         isz_i   = isz[is_no];
289d5b485abSSatish Balay         data_i  = data[is_no];
290d5b485abSSatish Balay         table_i = table[is_no];
291d5b485abSSatish Balay         for (k = 0; k < max; k++, ct1++) {
292d5b485abSSatish Balay           row = rbuf2_i[ct1];
29326fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row;
294d5b485abSSatish Balay         }
295d5b485abSSatish Balay         isz[is_no] = isz_i;
296d5b485abSSatish Balay       }
297d5b485abSSatish Balay     }
2989566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE));
299d5b485abSSatish Balay   }
300d5b485abSSatish Balay 
301d5b485abSSatish Balay   for (i = 0; i < imax; ++i) {
3029566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(iscomms[i], isz[i], data[i], PETSC_COPY_VALUES, is + i));
3039566063dSJacob Faibussowitsch     PetscCall(PetscCommDestroy(&iscomms[i]));
304d5b485abSSatish Balay   }
305d5b485abSSatish Balay 
3069566063dSJacob Faibussowitsch   PetscCall(PetscFree(iscomms));
3079566063dSJacob Faibussowitsch   PetscCall(PetscFree(onodes2));
3089566063dSJacob Faibussowitsch   PetscCall(PetscFree(olengths2));
309c7dd2924SSatish Balay 
3109566063dSJacob Faibussowitsch   PetscCall(PetscFree(pa));
311175c2bc5SJunchao Zhang   if (rbuf2) {
3129566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf2[0]));
3139566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf2));
314175c2bc5SJunchao Zhang   }
3159566063dSJacob Faibussowitsch   PetscCall(PetscFree(s_waits1));
3169566063dSJacob Faibussowitsch   PetscCall(PetscFree(r_waits1));
3179566063dSJacob Faibussowitsch   PetscCall(PetscFree(s_waits2));
3189566063dSJacob Faibussowitsch   PetscCall(PetscFree(r_waits2));
3199566063dSJacob Faibussowitsch   PetscCall(PetscFree5(table, data, isz, d_p, t_p));
320175c2bc5SJunchao Zhang   if (xdata) {
3219566063dSJacob Faibussowitsch     PetscCall(PetscFree(xdata[0]));
3229566063dSJacob Faibussowitsch     PetscCall(PetscFree(xdata));
323175c2bc5SJunchao Zhang   }
3249566063dSJacob Faibussowitsch   PetscCall(PetscFree(isz1));
3253a40ed3dSBarry Smith   PetscFunctionReturn(0);
326d5b485abSSatish Balay }
327d5b485abSSatish Balay 
328d5b485abSSatish Balay /*
329df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
330d5b485abSSatish Balay        the work on the local processor.
331d5b485abSSatish Balay 
332d5b485abSSatish Balay      Inputs:
333df36731bSBarry Smith       C      - MAT_MPIBAIJ;
334d5b485abSSatish Balay       imax - total no of index sets processed at a time;
335df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
336d5b485abSSatish Balay 
337d5b485abSSatish Balay      Output:
3380e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
339d5b485abSSatish Balay                to each index set;
340d5b485abSSatish Balay       data   - pointer to the solutions
341d5b485abSSatish Balay */
342d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C, PetscInt imax, PetscBT *table, PetscInt *isz, PetscInt **data)
343d71ae5a4SJacob Faibussowitsch {
344df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data;
345d5b485abSSatish Balay   Mat          A = c->A, B = c->B;
346df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
347b24ad042SBarry Smith   PetscInt     start, end, val, max, rstart, cstart, *ai, *aj;
348b24ad042SBarry Smith   PetscInt    *bi, *bj, *garray, i, j, k, row, *data_i, isz_i;
349f1af5d2fSBarry Smith   PetscBT      table_i;
350d5b485abSSatish Balay 
3513a40ed3dSBarry Smith   PetscFunctionBegin;
352899cda47SBarry Smith   rstart = c->rstartbs;
353899cda47SBarry Smith   cstart = c->cstartbs;
354d5b485abSSatish Balay   ai     = a->i;
355df36731bSBarry Smith   aj     = a->j;
356d5b485abSSatish Balay   bi     = b->i;
357df36731bSBarry Smith   bj     = b->j;
358d5b485abSSatish Balay   garray = c->garray;
359d5b485abSSatish Balay 
360d5b485abSSatish Balay   for (i = 0; i < imax; i++) {
361d5b485abSSatish Balay     data_i  = data[i];
362d5b485abSSatish Balay     table_i = table[i];
363d5b485abSSatish Balay     isz_i   = isz[i];
364d5b485abSSatish Balay     for (j = 0, max = isz[i]; j < max; j++) {
365d5b485abSSatish Balay       row   = data_i[j] - rstart;
366d5b485abSSatish Balay       start = ai[row];
367d5b485abSSatish Balay       end   = ai[row + 1];
368d5b485abSSatish Balay       for (k = start; k < end; k++) { /* Amat */
369df36731bSBarry Smith         val = aj[k] + cstart;
37026fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val;
371d5b485abSSatish Balay       }
372d5b485abSSatish Balay       start = bi[row];
373d5b485abSSatish Balay       end   = bi[row + 1];
374d5b485abSSatish Balay       for (k = start; k < end; k++) { /* Bmat */
375df36731bSBarry Smith         val = garray[bj[k]];
37626fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val;
377d5b485abSSatish Balay       }
378d5b485abSSatish Balay     }
379d5b485abSSatish Balay     isz[i] = isz_i;
380d5b485abSSatish Balay   }
3813a40ed3dSBarry Smith   PetscFunctionReturn(0);
382d5b485abSSatish Balay }
383d5b485abSSatish Balay /*
3842d4ee042Sprj-       MatIncreaseOverlap_MPIBAIJ_Receive - Process the received messages,
385d5b485abSSatish Balay          and return the output
386d5b485abSSatish Balay 
387d5b485abSSatish Balay          Input:
388d5b485abSSatish Balay            C    - the matrix
389d5b485abSSatish Balay            nrqr - no of messages being processed.
3902d4ee042Sprj-            rbuf - an array of pointers to the received requests
391d5b485abSSatish Balay 
392d5b485abSSatish Balay          Output:
393d5b485abSSatish Balay            xdata - array of messages to be sent back
394d5b485abSSatish Balay            isz1  - size of each message
395d5b485abSSatish Balay 
396a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
397d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
398a5b23f4aSJose E. Roman rather than all previous rows as it is now where a single large chunk of
399d5b485abSSatish Balay memory is used.
400d5b485abSSatish Balay 
401d5b485abSSatish Balay */
402d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C, PetscInt nrqr, PetscInt **rbuf, PetscInt **xdata, PetscInt *isz1)
403d71ae5a4SJacob Faibussowitsch {
404df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data;
405d5b485abSSatish Balay   Mat          A = c->A, B = c->B;
406df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
407b24ad042SBarry Smith   PetscInt     rstart, cstart, *ai, *aj, *bi, *bj, *garray, i, j, k;
408b24ad042SBarry Smith   PetscInt     row, total_sz, ct, ct1, ct2, ct3, mem_estimate, oct2, l, start, end;
409d2a212eaSBarry Smith   PetscInt     val, max1, max2, Mbs, no_malloc = 0, *tmp, new_estimate, ctr;
410b24ad042SBarry Smith   PetscInt    *rbuf_i, kmax, rbuf_0;
411f1af5d2fSBarry Smith   PetscBT      xtable;
412d5b485abSSatish Balay 
4133a40ed3dSBarry Smith   PetscFunctionBegin;
414df36731bSBarry Smith   Mbs    = c->Mbs;
415899cda47SBarry Smith   rstart = c->rstartbs;
416899cda47SBarry Smith   cstart = c->cstartbs;
417d5b485abSSatish Balay   ai     = a->i;
418df36731bSBarry Smith   aj     = a->j;
419d5b485abSSatish Balay   bi     = b->i;
420df36731bSBarry Smith   bj     = b->j;
421d5b485abSSatish Balay   garray = c->garray;
422d5b485abSSatish Balay 
423d5b485abSSatish Balay   for (i = 0, ct = 0, total_sz = 0; i < nrqr; ++i) {
424d5b485abSSatish Balay     rbuf_i = rbuf[i];
425d5b485abSSatish Balay     rbuf_0 = rbuf_i[0];
426d5b485abSSatish Balay     ct += rbuf_0;
42726fbe8dcSKarl Rupp     for (j = 1; j <= rbuf_0; j++) total_sz += rbuf_i[2 * j];
428d5b485abSSatish Balay   }
429d5b485abSSatish Balay 
430701b8814SBarry Smith   if (c->Mbs) max1 = ct * (a->nz + b->nz) / c->Mbs;
431701b8814SBarry Smith   else max1 = 1;
432d5b485abSSatish Balay   mem_estimate = 3 * ((total_sz > max1 ? total_sz : max1) + 1);
433175c2bc5SJunchao Zhang   if (nrqr) {
4349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mem_estimate, &xdata[0]));
435d5b485abSSatish Balay     ++no_malloc;
436175c2bc5SJunchao Zhang   }
4379566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(Mbs, &xtable));
4389566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(isz1, nrqr));
439d5b485abSSatish Balay 
440d5b485abSSatish Balay   ct3 = 0;
441d5b485abSSatish Balay   for (i = 0; i < nrqr; i++) { /* for easch mesg from proc i */
442d5b485abSSatish Balay     rbuf_i = rbuf[i];
443d5b485abSSatish Balay     rbuf_0 = rbuf_i[0];
444d5b485abSSatish Balay     ct1    = 2 * rbuf_0 + 1;
445d5b485abSSatish Balay     ct2    = ct1;
446d5b485abSSatish Balay     ct3 += ct1;
447d5b485abSSatish Balay     for (j = 1; j <= rbuf_0; j++) { /* for each IS from proc i*/
4489566063dSJacob Faibussowitsch       PetscCall(PetscBTMemzero(Mbs, xtable));
449d5b485abSSatish Balay       oct2 = ct2;
450d5b485abSSatish Balay       kmax = rbuf_i[2 * j];
451d5b485abSSatish Balay       for (k = 0; k < kmax; k++, ct1++) {
452d5b485abSSatish Balay         row = rbuf_i[ct1];
453f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable, row)) {
454d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
455b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5 * mem_estimate) + 1;
4569566063dSJacob Faibussowitsch             PetscCall(PetscMalloc1(new_estimate, &tmp));
4579566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate));
4589566063dSJacob Faibussowitsch             PetscCall(PetscFree(xdata[0]));
459d5b485abSSatish Balay             xdata[0]     = tmp;
4609371c9d4SSatish Balay             mem_estimate = new_estimate;
4619371c9d4SSatish Balay             ++no_malloc;
46226fbe8dcSKarl Rupp             for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1];
463d5b485abSSatish Balay           }
464d5b485abSSatish Balay           xdata[i][ct2++] = row;
465d5b485abSSatish Balay           ct3++;
466d5b485abSSatish Balay         }
467d5b485abSSatish Balay       }
468d5b485abSSatish Balay       for (k = oct2, max2 = ct2; k < max2; k++) {
469d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
470d5b485abSSatish Balay         start = ai[row];
471d5b485abSSatish Balay         end   = ai[row + 1];
472d5b485abSSatish Balay         for (l = start; l < end; l++) {
473df36731bSBarry Smith           val = aj[l] + cstart;
474f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable, val)) {
475d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
476b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5 * mem_estimate) + 1;
4779566063dSJacob Faibussowitsch               PetscCall(PetscMalloc1(new_estimate, &tmp));
4789566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate));
4799566063dSJacob Faibussowitsch               PetscCall(PetscFree(xdata[0]));
480d5b485abSSatish Balay               xdata[0]     = tmp;
4819371c9d4SSatish Balay               mem_estimate = new_estimate;
4829371c9d4SSatish Balay               ++no_malloc;
48326fbe8dcSKarl Rupp               for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1];
484d5b485abSSatish Balay             }
485d5b485abSSatish Balay             xdata[i][ct2++] = val;
486d5b485abSSatish Balay             ct3++;
487d5b485abSSatish Balay           }
488d5b485abSSatish Balay         }
489d5b485abSSatish Balay         start = bi[row];
490d5b485abSSatish Balay         end   = bi[row + 1];
491d5b485abSSatish Balay         for (l = start; l < end; l++) {
492df36731bSBarry Smith           val = garray[bj[l]];
493f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable, val)) {
494d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
495b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5 * mem_estimate) + 1;
4969566063dSJacob Faibussowitsch               PetscCall(PetscMalloc1(new_estimate, &tmp));
4979566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate));
4989566063dSJacob Faibussowitsch               PetscCall(PetscFree(xdata[0]));
499d5b485abSSatish Balay               xdata[0]     = tmp;
5009371c9d4SSatish Balay               mem_estimate = new_estimate;
5019371c9d4SSatish Balay               ++no_malloc;
50226fbe8dcSKarl Rupp               for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1];
503d5b485abSSatish Balay             }
504d5b485abSSatish Balay             xdata[i][ct2++] = val;
505d5b485abSSatish Balay             ct3++;
506d5b485abSSatish Balay           }
507d5b485abSSatish Balay         }
508d5b485abSSatish Balay       }
509d5b485abSSatish Balay       /* Update the header*/
510d5b485abSSatish Balay       xdata[i][2 * j]     = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
511d5b485abSSatish Balay       xdata[i][2 * j - 1] = rbuf_i[2 * j - 1];
512d5b485abSSatish Balay     }
513d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
514175c2bc5SJunchao Zhang     if (i + 1 < nrqr) xdata[i + 1] = xdata[i] + ct2;
515d5b485abSSatish Balay     isz1[i] = ct2; /* size of each message */
516d5b485abSSatish Balay   }
5179566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&xtable));
5189566063dSJacob Faibussowitsch   PetscCall(PetscInfo(C, "Allocated %" PetscInt_FMT " bytes, required %" PetscInt_FMT ", no of mallocs = %" PetscInt_FMT "\n", mem_estimate, ct3, no_malloc));
5193a40ed3dSBarry Smith   PetscFunctionReturn(0);
520d5b485abSSatish Balay }
521d5b485abSSatish Balay 
522d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submat[])
523d71ae5a4SJacob Faibussowitsch {
52417df9f7cSHong Zhang   IS          *isrow_block, *iscol_block;
525cf4f527aSSatish Balay   Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data;
526121971b2SHong Zhang   PetscInt     nmax, nstages_local, nstages, i, pos, max_no, N = C->cmap->N, bs = C->rmap->bs;
5275dc98d6aSHong Zhang   Mat_SeqBAIJ *subc;
5285c39f6d9SHong Zhang   Mat_SubSppt *smat;
529a2feddc5SSatish Balay 
5303a40ed3dSBarry Smith   PetscFunctionBegin;
53129dcf524SDmitry Karpeev   /* The compression and expansion should be avoided. Doesn't point
53229dcf524SDmitry Karpeev      out errors, might change the indices, hence buggey */
5339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(ismax + 1, &isrow_block, ismax + 1, &iscol_block));
5349566063dSJacob Faibussowitsch   PetscCall(ISCompressIndicesGeneral(N, C->rmap->n, bs, ismax, isrow, isrow_block));
5359566063dSJacob Faibussowitsch   PetscCall(ISCompressIndicesGeneral(N, C->cmap->n, bs, ismax, iscol, iscol_block));
536cf4f527aSSatish Balay 
537cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
5384698041cSHong Zhang   if (!C->cmap->N) nmax = 20 * 1000000 / sizeof(PetscInt);
5394698041cSHong Zhang   else nmax = 20 * 1000000 / (c->Nbs * sizeof(PetscInt));
540329f5518SBarry Smith   if (!nmax) nmax = 1;
5415dc98d6aSHong Zhang 
5425dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
5434698041cSHong Zhang     nstages_local = ismax / nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */
544cf4f527aSSatish Balay 
545653e4784SBarry Smith     /* Make sure every processor loops through the nstages */
5461c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&nstages_local, &nstages, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)C)));
5475dc98d6aSHong Zhang 
5484698041cSHong Zhang     /* Allocate memory to hold all the submatrices and dummy submatrices */
5499566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(ismax + nstages, submat));
5504698041cSHong Zhang   } else { /* MAT_REUSE_MATRIX */
5514698041cSHong Zhang     if (ismax) {
5525dc98d6aSHong Zhang       subc = (Mat_SeqBAIJ *)((*submat)[0]->data);
5535dc98d6aSHong Zhang       smat = subc->submatis1;
5544698041cSHong Zhang     } else { /* (*submat)[0] is a dummy matrix */
5555c39f6d9SHong Zhang       smat = (Mat_SubSppt *)(*submat)[0]->data;
5564698041cSHong Zhang     }
55728b400f6SJacob Faibussowitsch     PetscCheck(smat, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat");
5585dc98d6aSHong Zhang     nstages = smat->nstages;
5595dc98d6aSHong Zhang   }
5605dc98d6aSHong Zhang 
561cf4f527aSSatish Balay   for (i = 0, pos = 0; i < nstages; i++) {
562cf4f527aSSatish Balay     if (pos + nmax <= ismax) max_no = nmax;
5634e195584SJunchao Zhang     else if (pos >= ismax) max_no = 0;
564cf4f527aSSatish Balay     else max_no = ismax - pos;
565a42ab0d6SHong Zhang 
5669566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices_MPIBAIJ_local(C, max_no, isrow_block + pos, iscol_block + pos, scall, *submat + pos));
5674e195584SJunchao Zhang     if (!max_no) {
5684e195584SJunchao Zhang       if (scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */
5695c39f6d9SHong Zhang         smat          = (Mat_SubSppt *)(*submat)[pos]->data;
5704698041cSHong Zhang         smat->nstages = nstages;
5714698041cSHong Zhang       }
5724e195584SJunchao Zhang       pos++; /* advance to next dummy matrix if any */
5734e195584SJunchao Zhang     } else pos += max_no;
574cf4f527aSSatish Balay   }
57536f4e84dSSatish Balay 
5765dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX && ismax) {
5775dc98d6aSHong Zhang     /* save nstages for reuse */
5785dc98d6aSHong Zhang     subc          = (Mat_SeqBAIJ *)((*submat)[0]->data);
5795dc98d6aSHong Zhang     smat          = subc->submatis1;
5805dc98d6aSHong Zhang     smat->nstages = nstages;
5815dc98d6aSHong Zhang   }
5825dc98d6aSHong Zhang 
58336f4e84dSSatish Balay   for (i = 0; i < ismax; i++) {
5849566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isrow_block[i]));
5859566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iscol_block[i]));
58636f4e84dSSatish Balay   }
5879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(isrow_block, iscol_block));
5883a40ed3dSBarry Smith   PetscFunctionReturn(0);
589a2feddc5SSatish Balay }
590a2feddc5SSatish Balay 
591233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
592d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
593d71ae5a4SJacob Faibussowitsch {
594e005ede5SBarry Smith   PetscInt    nGlobalNd = proc_gnode[size];
5954dc2109aSBarry Smith   PetscMPIInt fproc;
596233c2ff6SSatish Balay 
597233c2ff6SSatish Balay   PetscFunctionBegin;
5989566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)), &fproc));
59923ce1328SBarry Smith   if (fproc > size) fproc = size;
600e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc + 1]) {
601e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
602233c2ff6SSatish Balay     else fproc++;
603233c2ff6SSatish Balay   }
604e005ede5SBarry Smith   *rank = fproc;
605233c2ff6SSatish Balay   PetscFunctionReturn(0);
606233c2ff6SSatish Balay }
607233c2ff6SSatish Balay #endif
608233c2ff6SSatish Balay 
609a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
610b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
611d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submats)
612d71ae5a4SJacob Faibussowitsch {
61317df9f7cSHong Zhang   Mat_MPIBAIJ     *c = (Mat_MPIBAIJ *)C->data;
61417df9f7cSHong Zhang   Mat              A = c->A;
61517df9f7cSHong Zhang   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)c->B->data, *subc;
61617df9f7cSHong Zhang   const PetscInt **icol, **irow;
61717df9f7cSHong Zhang   PetscInt        *nrow, *ncol, start;
61817df9f7cSHong Zhang   PetscMPIInt      rank, size, tag0, tag2, tag3, tag4, *w1, *w2, *w3, *w4, nrqr;
619ddb3d6bcSHong Zhang   PetscInt       **sbuf1, **sbuf2, *sbuf2_i, i, j, k, l, ct1, ct2, **rbuf1, row, proc = -1;
62017df9f7cSHong Zhang   PetscInt         nrqs = 0, msz, **ptr = NULL, *req_size = NULL, *ctr = NULL, *pa, *tmp = NULL, tcol;
62117df9f7cSHong Zhang   PetscInt       **rbuf3 = NULL, *req_source1 = NULL, *req_source2, **sbuf_aj, **rbuf2 = NULL, max1, max2;
62217df9f7cSHong Zhang   PetscInt       **lens, is_no, ncols, *cols, mat_i, *mat_j, tmp2, jmax;
62317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
624*eec179cfSJacob Faibussowitsch   PetscHMapI *cmap, cmap_i = NULL, *rmap, rmap_i;
62517df9f7cSHong Zhang #else
62617df9f7cSHong Zhang   PetscInt **cmap, *cmap_i = NULL, **rmap, *rmap_i;
62717df9f7cSHong Zhang #endif
62896cff833SHong Zhang   const PetscInt *irow_i, *icol_i;
62917df9f7cSHong Zhang   PetscInt        ctr_j, *sbuf1_j, *sbuf_aj_i, *rbuf1_i, kmax, *lens_i;
63017df9f7cSHong Zhang   MPI_Request    *s_waits1, *r_waits1, *s_waits2, *r_waits2, *r_waits3;
63117df9f7cSHong Zhang   MPI_Request    *r_waits4, *s_waits3, *s_waits4;
63217df9f7cSHong Zhang   MPI_Comm        comm;
63364f5bb2dSSatish Balay   PetscScalar   **rbuf4, *rbuf4_i = NULL, **sbuf_aa, *vals, *mat_a = NULL, *imat_a = NULL, *sbuf_aa_i;
63417df9f7cSHong Zhang   PetscMPIInt    *onodes1, *olengths1, end;
63580d31651SHong Zhang   PetscInt      **row2proc, *row2proc_i, *imat_ilen, *imat_j, *imat_i;
6365c39f6d9SHong Zhang   Mat_SubSppt    *smat_i;
63717df9f7cSHong Zhang   PetscBool      *issorted, colflag, iscsorted = PETSC_TRUE;
63817df9f7cSHong Zhang   PetscInt       *sbuf1_i, *rbuf2_i, *rbuf3_i, ilen;
639a42ab0d6SHong Zhang   PetscInt        bs = C->rmap->bs, bs2 = c->bs2, rstart = c->rstartbs;
640a42ab0d6SHong Zhang   PetscBool       ijonly = c->ijonly; /* private flag indicates only matrix data structures are requested */
641a42ab0d6SHong Zhang   PetscInt        nzA, nzB, *a_i = a->i, *b_i = b->i, *a_j = a->j, *b_j = b->j, ctmp, imark, *cworkA, *cworkB;
642afb3cbd8SHong Zhang   PetscScalar    *vworkA = NULL, *vworkB = NULL, *a_a = a->a, *b_a = b->a;
643a42ab0d6SHong Zhang   PetscInt        cstart = c->cstartbs, *bmap = c->garray;
64480d31651SHong Zhang   PetscBool      *allrows, *allcolumns;
645a42ab0d6SHong Zhang 
64617df9f7cSHong Zhang   PetscFunctionBegin;
6479566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
64817df9f7cSHong Zhang   size = c->size;
64917df9f7cSHong Zhang   rank = c->rank;
65017df9f7cSHong Zhang 
6519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(ismax, &row2proc, ismax, &cmap, ismax, &rmap, ismax + 1, &allcolumns, ismax, &allrows));
6529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(ismax, (PetscInt ***)&irow, ismax, (PetscInt ***)&icol, ismax, &nrow, ismax, &ncol, ismax, &issorted));
65317df9f7cSHong Zhang 
65417df9f7cSHong Zhang   for (i = 0; i < ismax; i++) {
6559566063dSJacob Faibussowitsch     PetscCall(ISSorted(iscol[i], &issorted[i]));
65680d31651SHong Zhang     if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */
6579566063dSJacob Faibussowitsch     PetscCall(ISSorted(isrow[i], &issorted[i]));
65817df9f7cSHong Zhang 
659ec3c739cSHong Zhang     /* Check for special case: allcolumns */
6609566063dSJacob Faibussowitsch     PetscCall(ISIdentity(iscol[i], &colflag));
6619566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(iscol[i], &ncol[i]));
6625dd0c08cSHong Zhang 
6631abc2f7bSHong Zhang     if (colflag && ncol[i] == c->Nbs) {
6641abc2f7bSHong Zhang       allcolumns[i] = PETSC_TRUE;
6651abc2f7bSHong Zhang       icol[i]       = NULL;
6661abc2f7bSHong Zhang     } else {
66717df9f7cSHong Zhang       allcolumns[i] = PETSC_FALSE;
6689566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(iscol[i], &icol[i]));
6691abc2f7bSHong Zhang     }
670ec3c739cSHong Zhang 
671ec3c739cSHong Zhang     /* Check for special case: allrows */
6729566063dSJacob Faibussowitsch     PetscCall(ISIdentity(isrow[i], &colflag));
6739566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(isrow[i], &nrow[i]));
674ec3c739cSHong Zhang     if (colflag && nrow[i] == c->Mbs) {
675ec3c739cSHong Zhang       allrows[i] = PETSC_TRUE;
676ec3c739cSHong Zhang       irow[i]    = NULL;
677ec3c739cSHong Zhang     } else {
678ec3c739cSHong Zhang       allrows[i] = PETSC_FALSE;
6799566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(isrow[i], &irow[i]));
680ec3c739cSHong Zhang     }
68117df9f7cSHong Zhang   }
68217df9f7cSHong Zhang 
68317df9f7cSHong Zhang   if (scall == MAT_REUSE_MATRIX) {
68417df9f7cSHong Zhang     /* Assumes new rows are same length as the old rows */
68517df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
68617df9f7cSHong Zhang       subc = (Mat_SeqBAIJ *)(submats[i]->data);
687aed4548fSBarry Smith       PetscCheck(subc->mbs == nrow[i] && subc->nbs == ncol[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size");
68817df9f7cSHong Zhang 
68917df9f7cSHong Zhang       /* Initial matrix as if empty */
6909566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(subc->ilen, subc->mbs));
69117df9f7cSHong Zhang 
69217df9f7cSHong Zhang       /* Initial matrix as if empty */
69317df9f7cSHong Zhang       submats[i]->factortype = C->factortype;
69417df9f7cSHong Zhang 
69517df9f7cSHong Zhang       smat_i = subc->submatis1;
69617df9f7cSHong Zhang 
69717df9f7cSHong Zhang       nrqs        = smat_i->nrqs;
69817df9f7cSHong Zhang       nrqr        = smat_i->nrqr;
69917df9f7cSHong Zhang       rbuf1       = smat_i->rbuf1;
70017df9f7cSHong Zhang       rbuf2       = smat_i->rbuf2;
70117df9f7cSHong Zhang       rbuf3       = smat_i->rbuf3;
70217df9f7cSHong Zhang       req_source2 = smat_i->req_source2;
70317df9f7cSHong Zhang 
70417df9f7cSHong Zhang       sbuf1 = smat_i->sbuf1;
70517df9f7cSHong Zhang       sbuf2 = smat_i->sbuf2;
70617df9f7cSHong Zhang       ptr   = smat_i->ptr;
70717df9f7cSHong Zhang       tmp   = smat_i->tmp;
70817df9f7cSHong Zhang       ctr   = smat_i->ctr;
70917df9f7cSHong Zhang 
71017df9f7cSHong Zhang       pa          = smat_i->pa;
71117df9f7cSHong Zhang       req_size    = smat_i->req_size;
71217df9f7cSHong Zhang       req_source1 = smat_i->req_source1;
71317df9f7cSHong Zhang 
71417df9f7cSHong Zhang       allcolumns[i] = smat_i->allcolumns;
715ec3c739cSHong Zhang       allrows[i]    = smat_i->allrows;
71617df9f7cSHong Zhang       row2proc[i]   = smat_i->row2proc;
71717df9f7cSHong Zhang       rmap[i]       = smat_i->rmap;
71817df9f7cSHong Zhang       cmap[i]       = smat_i->cmap;
71917df9f7cSHong Zhang     }
7204698041cSHong Zhang 
7214698041cSHong Zhang     if (!ismax) { /* Get dummy submatrices and retrieve struct submatis1 */
72208401ef6SPierre Jolivet       PetscCheck(submats[0], PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "submats are null, cannot reuse");
7235c39f6d9SHong Zhang       smat_i = (Mat_SubSppt *)submats[0]->data;
7244698041cSHong Zhang 
7254698041cSHong Zhang       nrqs        = smat_i->nrqs;
7264698041cSHong Zhang       nrqr        = smat_i->nrqr;
7274698041cSHong Zhang       rbuf1       = smat_i->rbuf1;
7284698041cSHong Zhang       rbuf2       = smat_i->rbuf2;
7294698041cSHong Zhang       rbuf3       = smat_i->rbuf3;
7304698041cSHong Zhang       req_source2 = smat_i->req_source2;
7314698041cSHong Zhang 
7324698041cSHong Zhang       sbuf1 = smat_i->sbuf1;
7334698041cSHong Zhang       sbuf2 = smat_i->sbuf2;
7344698041cSHong Zhang       ptr   = smat_i->ptr;
7354698041cSHong Zhang       tmp   = smat_i->tmp;
7364698041cSHong Zhang       ctr   = smat_i->ctr;
7374698041cSHong Zhang 
7384698041cSHong Zhang       pa          = smat_i->pa;
7394698041cSHong Zhang       req_size    = smat_i->req_size;
7404698041cSHong Zhang       req_source1 = smat_i->req_source1;
7414698041cSHong Zhang 
742c67fe082SHong Zhang       allcolumns[0] = PETSC_FALSE;
7434698041cSHong Zhang     }
74417df9f7cSHong Zhang   } else { /* scall == MAT_INITIAL_MATRIX */
74517df9f7cSHong Zhang     /* Get some new tags to keep the communication clean */
7469566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2));
7479566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag3));
74817df9f7cSHong Zhang 
74917df9f7cSHong Zhang     /* evaluate communication - mesg to who, length of mesg, and buffer space
75017df9f7cSHong Zhang      required. Based on this, buffers are allocated, and data copied into them*/
7519566063dSJacob Faibussowitsch     PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4)); /* mesg size, initialize work vectors */
75217df9f7cSHong Zhang 
75317df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
75417df9f7cSHong Zhang       jmax   = nrow[i];
75517df9f7cSHong Zhang       irow_i = irow[i];
75617df9f7cSHong Zhang 
7579566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jmax, &row2proc_i));
75817df9f7cSHong Zhang       row2proc[i] = row2proc_i;
75917df9f7cSHong Zhang 
76017df9f7cSHong Zhang       if (issorted[i]) proc = 0;
76117df9f7cSHong Zhang       for (j = 0; j < jmax; j++) {
76217df9f7cSHong Zhang         if (!issorted[i]) proc = 0;
763ec3c739cSHong Zhang         if (allrows[i]) row = j;
764ec3c739cSHong Zhang         else row = irow_i[j];
765ec3c739cSHong Zhang 
766a42ab0d6SHong Zhang         while (row >= c->rangebs[proc + 1]) proc++;
76717df9f7cSHong Zhang         w4[proc]++;
76817df9f7cSHong Zhang         row2proc_i[j] = proc; /* map row index to proc */
76917df9f7cSHong Zhang       }
77017df9f7cSHong Zhang       for (j = 0; j < size; j++) {
7719371c9d4SSatish Balay         if (w4[j]) {
7729371c9d4SSatish Balay           w1[j] += w4[j];
7739371c9d4SSatish Balay           w3[j]++;
7749371c9d4SSatish Balay           w4[j] = 0;
7759371c9d4SSatish Balay         }
77617df9f7cSHong Zhang       }
77717df9f7cSHong Zhang     }
77817df9f7cSHong Zhang 
77917df9f7cSHong Zhang     nrqs     = 0; /* no of outgoing messages */
78017df9f7cSHong Zhang     msz      = 0; /* total mesg length (for all procs) */
78117df9f7cSHong Zhang     w1[rank] = 0; /* no mesg sent to self */
78217df9f7cSHong Zhang     w3[rank] = 0;
78317df9f7cSHong Zhang     for (i = 0; i < size; i++) {
7849371c9d4SSatish Balay       if (w1[i]) {
7859371c9d4SSatish Balay         w2[i] = 1;
7869371c9d4SSatish Balay         nrqs++;
7879371c9d4SSatish Balay       } /* there exists a message to proc i */
78817df9f7cSHong Zhang     }
7899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &pa)); /*(proc -array)*/
79017df9f7cSHong Zhang     for (i = 0, j = 0; i < size; i++) {
7919371c9d4SSatish Balay       if (w1[i]) {
7929371c9d4SSatish Balay         pa[j] = i;
7939371c9d4SSatish Balay         j++;
7949371c9d4SSatish Balay       }
79517df9f7cSHong Zhang     }
79617df9f7cSHong Zhang 
79717df9f7cSHong Zhang     /* Each message would have a header = 1 + 2*(no of IS) + data */
79817df9f7cSHong Zhang     for (i = 0; i < nrqs; i++) {
79917df9f7cSHong Zhang       j = pa[i];
80017df9f7cSHong Zhang       w1[j] += w2[j] + 2 * w3[j];
80117df9f7cSHong Zhang       msz += w1[j];
80217df9f7cSHong Zhang     }
8039566063dSJacob Faibussowitsch     PetscCall(PetscInfo(0, "Number of outgoing messages %" PetscInt_FMT " Total message length %" PetscInt_FMT "\n", nrqs, msz));
80417df9f7cSHong Zhang 
80517df9f7cSHong Zhang     /* Determine the number of messages to expect, their lengths, from from-ids */
8069566063dSJacob Faibussowitsch     PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr));
8079566063dSJacob Faibussowitsch     PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1));
80817df9f7cSHong Zhang 
80917df9f7cSHong Zhang     /* Now post the Irecvs corresponding to these messages */
8109566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag0));
8119566063dSJacob Faibussowitsch     PetscCall(PetscPostIrecvInt(comm, tag0, nrqr, onodes1, olengths1, &rbuf1, &r_waits1));
81217df9f7cSHong Zhang 
81317df9f7cSHong Zhang     /* Allocate Memory for outgoing messages */
8149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(size, &sbuf1, size, &ptr, 2 * msz, &tmp, size, &ctr));
8159566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(sbuf1, size));
8169566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ptr, size));
81717df9f7cSHong Zhang 
81817df9f7cSHong Zhang     {
81917df9f7cSHong Zhang       PetscInt *iptr = tmp;
82017df9f7cSHong Zhang       k              = 0;
82117df9f7cSHong Zhang       for (i = 0; i < nrqs; i++) {
82217df9f7cSHong Zhang         j = pa[i];
82317df9f7cSHong Zhang         iptr += k;
82417df9f7cSHong Zhang         sbuf1[j] = iptr;
82517df9f7cSHong Zhang         k        = w1[j];
82617df9f7cSHong Zhang       }
82717df9f7cSHong Zhang     }
82817df9f7cSHong Zhang 
82917df9f7cSHong Zhang     /* Form the outgoing messages. Initialize the header space */
83017df9f7cSHong Zhang     for (i = 0; i < nrqs; i++) {
83117df9f7cSHong Zhang       j           = pa[i];
83217df9f7cSHong Zhang       sbuf1[j][0] = 0;
8339566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(sbuf1[j] + 1, 2 * w3[j]));
83417df9f7cSHong Zhang       ptr[j] = sbuf1[j] + 2 * w3[j] + 1;
83517df9f7cSHong Zhang     }
83617df9f7cSHong Zhang 
83717df9f7cSHong Zhang     /* Parse the isrow and copy data into outbuf */
83817df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
83917df9f7cSHong Zhang       row2proc_i = row2proc[i];
8409566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctr, size));
84117df9f7cSHong Zhang       irow_i = irow[i];
84217df9f7cSHong Zhang       jmax   = nrow[i];
84317df9f7cSHong Zhang       for (j = 0; j < jmax; j++) { /* parse the indices of each IS */
84417df9f7cSHong Zhang         proc = row2proc_i[j];
845ec3c739cSHong Zhang         if (allrows[i]) row = j;
846ec3c739cSHong Zhang         else row = irow_i[j];
847ec3c739cSHong Zhang 
84817df9f7cSHong Zhang         if (proc != rank) { /* copy to the outgoing buf*/
84917df9f7cSHong Zhang           ctr[proc]++;
850ec3c739cSHong Zhang           *ptr[proc] = row;
85117df9f7cSHong Zhang           ptr[proc]++;
85217df9f7cSHong Zhang         }
85317df9f7cSHong Zhang       }
85417df9f7cSHong Zhang       /* Update the headers for the current IS */
85517df9f7cSHong Zhang       for (j = 0; j < size; j++) { /* Can Optimise this loop too */
85617df9f7cSHong Zhang         if ((ctr_j = ctr[j])) {
85717df9f7cSHong Zhang           sbuf1_j            = sbuf1[j];
85817df9f7cSHong Zhang           k                  = ++sbuf1_j[0];
85917df9f7cSHong Zhang           sbuf1_j[2 * k]     = ctr_j;
86017df9f7cSHong Zhang           sbuf1_j[2 * k - 1] = i;
86117df9f7cSHong Zhang         }
86217df9f7cSHong Zhang       }
86317df9f7cSHong Zhang     }
86417df9f7cSHong Zhang 
86517df9f7cSHong Zhang     /*  Now  post the sends */
8669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &s_waits1));
86717df9f7cSHong Zhang     for (i = 0; i < nrqs; ++i) {
86817df9f7cSHong Zhang       j = pa[i];
8699566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sbuf1[j], w1[j], MPIU_INT, j, tag0, comm, s_waits1 + i));
87017df9f7cSHong Zhang     }
87117df9f7cSHong Zhang 
87217df9f7cSHong Zhang     /* Post Receives to capture the buffer size */
8739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &r_waits2));
8749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(nrqs, &req_source2, nrqs, &rbuf2, nrqs, &rbuf3));
875175c2bc5SJunchao Zhang     if (nrqs) rbuf2[0] = tmp + msz;
876ad540459SPierre Jolivet     for (i = 1; i < nrqs; ++i) rbuf2[i] = rbuf2[i - 1] + w1[pa[i - 1]];
87717df9f7cSHong Zhang     for (i = 0; i < nrqs; ++i) {
87817df9f7cSHong Zhang       j = pa[i];
8799566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(rbuf2[i], w1[j], MPIU_INT, j, tag2, comm, r_waits2 + i));
88017df9f7cSHong Zhang     }
88117df9f7cSHong Zhang 
88217df9f7cSHong Zhang     /* Send to other procs the buf size they should allocate */
88317df9f7cSHong Zhang     /* Receive messages*/
8849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &s_waits2));
8859566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(nrqr, &sbuf2, nrqr, &req_size, nrqr, &req_source1));
88617df9f7cSHong Zhang 
8879566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE));
88817df9f7cSHong Zhang     for (i = 0; i < nrqr; ++i) {
88917df9f7cSHong Zhang       req_size[i] = 0;
89017df9f7cSHong Zhang       rbuf1_i     = rbuf1[i];
89117df9f7cSHong Zhang       start       = 2 * rbuf1_i[0] + 1;
892175c2bc5SJunchao Zhang       end         = olengths1[i];
8939566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(end, &sbuf2[i]));
89417df9f7cSHong Zhang       sbuf2_i = sbuf2[i];
89517df9f7cSHong Zhang       for (j = start; j < end; j++) {
896ddb3d6bcSHong Zhang         row        = rbuf1_i[j] - rstart;
897ddb3d6bcSHong Zhang         ncols      = a_i[row + 1] - a_i[row] + b_i[row + 1] - b_i[row];
89817df9f7cSHong Zhang         sbuf2_i[j] = ncols;
89917df9f7cSHong Zhang         req_size[i] += ncols;
90017df9f7cSHong Zhang       }
901175c2bc5SJunchao Zhang       req_source1[i] = onodes1[i];
90217df9f7cSHong Zhang       /* form the header */
90317df9f7cSHong Zhang       sbuf2_i[0] = req_size[i];
90417df9f7cSHong Zhang       for (j = 1; j < start; j++) sbuf2_i[j] = rbuf1_i[j];
90517df9f7cSHong Zhang 
9069566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sbuf2_i, end, MPIU_INT, req_source1[i], tag2, comm, s_waits2 + i));
90717df9f7cSHong Zhang     }
908ddb3d6bcSHong Zhang 
9099566063dSJacob Faibussowitsch     PetscCall(PetscFree(onodes1));
9109566063dSJacob Faibussowitsch     PetscCall(PetscFree(olengths1));
911175c2bc5SJunchao Zhang 
9129566063dSJacob Faibussowitsch     PetscCall(PetscFree(r_waits1));
9139566063dSJacob Faibussowitsch     PetscCall(PetscFree4(w1, w2, w3, w4));
91417df9f7cSHong Zhang 
91517df9f7cSHong Zhang     /* Receive messages*/
9169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &r_waits3));
91717df9f7cSHong Zhang 
9189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqs, r_waits2, MPI_STATUSES_IGNORE));
91917df9f7cSHong Zhang     for (i = 0; i < nrqs; ++i) {
9209566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(rbuf2[i][0], &rbuf3[i]));
921175c2bc5SJunchao Zhang       req_source2[i] = pa[i];
9229566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(rbuf3[i], rbuf2[i][0], MPIU_INT, req_source2[i], tag3, comm, r_waits3 + i));
92317df9f7cSHong Zhang     }
9249566063dSJacob Faibussowitsch     PetscCall(PetscFree(r_waits2));
92517df9f7cSHong Zhang 
92617df9f7cSHong Zhang     /* Wait on sends1 and sends2 */
9279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE));
9289566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE));
9299566063dSJacob Faibussowitsch     PetscCall(PetscFree(s_waits1));
9309566063dSJacob Faibussowitsch     PetscCall(PetscFree(s_waits2));
93117df9f7cSHong Zhang 
93217df9f7cSHong Zhang     /* Now allocate sending buffers for a->j, and send them off */
9339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &sbuf_aj));
93417df9f7cSHong Zhang     for (i = 0, j = 0; i < nrqr; i++) j += req_size[i];
9359566063dSJacob Faibussowitsch     if (nrqr) PetscCall(PetscMalloc1(j, &sbuf_aj[0]));
93617df9f7cSHong Zhang     for (i = 1; i < nrqr; i++) sbuf_aj[i] = sbuf_aj[i - 1] + req_size[i - 1];
93717df9f7cSHong Zhang 
9389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &s_waits3));
93917df9f7cSHong Zhang     {
94017df9f7cSHong Zhang       for (i = 0; i < nrqr; i++) {
94117df9f7cSHong Zhang         rbuf1_i   = rbuf1[i];
94217df9f7cSHong Zhang         sbuf_aj_i = sbuf_aj[i];
94317df9f7cSHong Zhang         ct1       = 2 * rbuf1_i[0] + 1;
94417df9f7cSHong Zhang         ct2       = 0;
94517df9f7cSHong Zhang         for (j = 1, max1 = rbuf1_i[0]; j <= max1; j++) {
94617df9f7cSHong Zhang           kmax = rbuf1[i][2 * j];
94717df9f7cSHong Zhang           for (k = 0; k < kmax; k++, ct1++) {
94817df9f7cSHong Zhang             row    = rbuf1_i[ct1] - rstart;
9499371c9d4SSatish Balay             nzA    = a_i[row + 1] - a_i[row];
9509371c9d4SSatish Balay             nzB    = b_i[row + 1] - b_i[row];
95117df9f7cSHong Zhang             ncols  = nzA + nzB;
9529371c9d4SSatish Balay             cworkA = a_j + a_i[row];
9539371c9d4SSatish Balay             cworkB = b_j + b_i[row];
95417df9f7cSHong Zhang 
95517df9f7cSHong Zhang             /* load the column indices for this row into cols */
95617df9f7cSHong Zhang             cols = sbuf_aj_i + ct2;
95717df9f7cSHong Zhang             for (l = 0; l < nzB; l++) {
958a42ab0d6SHong Zhang               if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
959a42ab0d6SHong Zhang               else break;
96017df9f7cSHong Zhang             }
961a42ab0d6SHong Zhang             imark = l;
962ad540459SPierre Jolivet             for (l = 0; l < nzA; l++) cols[imark + l] = cstart + cworkA[l];
963a42ab0d6SHong Zhang             for (l = imark; l < nzB; l++) cols[nzA + l] = bmap[cworkB[l]];
96417df9f7cSHong Zhang             ct2 += ncols;
96517df9f7cSHong Zhang           }
96617df9f7cSHong Zhang         }
9679566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Isend(sbuf_aj_i, req_size[i], MPIU_INT, req_source1[i], tag3, comm, s_waits3 + i));
96817df9f7cSHong Zhang       }
96917df9f7cSHong Zhang     }
97017df9f7cSHong Zhang 
97117df9f7cSHong Zhang     /* create col map: global col of C -> local col of submatrices */
97217df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
97317df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
97417df9f7cSHong Zhang       if (!allcolumns[i]) {
975*eec179cfSJacob Faibussowitsch         PetscCall(PetscHMapICreateWithSize(ncol[i], cmap + i));
97617df9f7cSHong Zhang 
97717df9f7cSHong Zhang         jmax   = ncol[i];
97817df9f7cSHong Zhang         icol_i = icol[i];
97917df9f7cSHong Zhang         cmap_i = cmap[i];
980*eec179cfSJacob Faibussowitsch         for (j = 0; j < jmax; j++) PetscCall(PetscHMapISetWithMode(cmap[i], icol_i[j] + 1, j + 1, INSERT_VALUES));
98117df9f7cSHong Zhang       } else cmap[i] = NULL;
98217df9f7cSHong Zhang     }
98317df9f7cSHong Zhang #else
98417df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
98517df9f7cSHong Zhang       if (!allcolumns[i]) {
9869566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(c->Nbs, &cmap[i]));
98717df9f7cSHong Zhang         jmax   = ncol[i];
98817df9f7cSHong Zhang         icol_i = icol[i];
98917df9f7cSHong Zhang         cmap_i = cmap[i];
990a42ab0d6SHong Zhang         for (j = 0; j < jmax; j++) cmap_i[icol_i[j]] = j + 1;
99117df9f7cSHong Zhang       } else cmap[i] = NULL;
99217df9f7cSHong Zhang     }
99317df9f7cSHong Zhang #endif
99417df9f7cSHong Zhang 
99517df9f7cSHong Zhang     /* Create lens which is required for MatCreate... */
99617df9f7cSHong Zhang     for (i = 0, j = 0; i < ismax; i++) j += nrow[i];
9979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ismax, &lens));
99817df9f7cSHong Zhang 
99948a46eb9SPierre Jolivet     if (ismax) PetscCall(PetscCalloc1(j, &lens[0]));
100017df9f7cSHong Zhang     for (i = 1; i < ismax; i++) lens[i] = lens[i - 1] + nrow[i - 1];
100117df9f7cSHong Zhang 
100217df9f7cSHong Zhang     /* Update lens from local data */
100317df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
100417df9f7cSHong Zhang       row2proc_i = row2proc[i];
100517df9f7cSHong Zhang       jmax       = nrow[i];
100617df9f7cSHong Zhang       if (!allcolumns[i]) cmap_i = cmap[i];
100717df9f7cSHong Zhang       irow_i = irow[i];
100817df9f7cSHong Zhang       lens_i = lens[i];
100917df9f7cSHong Zhang       for (j = 0; j < jmax; j++) {
1010ec3c739cSHong Zhang         if (allrows[i]) row = j;
1011ec3c739cSHong Zhang         else row = irow_i[j]; /* global blocked row of C */
1012ec3c739cSHong Zhang 
101317df9f7cSHong Zhang         proc = row2proc_i[j];
101417df9f7cSHong Zhang         if (proc == rank) {
1015a42ab0d6SHong Zhang           /* Get indices from matA and then from matB */
101617df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1017a42ab0d6SHong Zhang           PetscInt tt;
1018a42ab0d6SHong Zhang #endif
1019a42ab0d6SHong Zhang           row    = row - rstart;
1020a42ab0d6SHong Zhang           nzA    = a_i[row + 1] - a_i[row];
1021a42ab0d6SHong Zhang           nzB    = b_i[row + 1] - b_i[row];
1022a42ab0d6SHong Zhang           cworkA = a_j + a_i[row];
1023a42ab0d6SHong Zhang           cworkB = b_j + b_i[row];
1024a42ab0d6SHong Zhang 
1025a42ab0d6SHong Zhang           if (!allcolumns[i]) {
1026a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1027a42ab0d6SHong Zhang             for (k = 0; k < nzA; k++) {
1028*eec179cfSJacob Faibussowitsch               PetscCall(PetscHMapIGetWithDefault(cmap_i, cstart + cworkA[k] + 1, 0, &tt));
1029a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1030a42ab0d6SHong Zhang             }
1031a42ab0d6SHong Zhang             for (k = 0; k < nzB; k++) {
1032*eec179cfSJacob Faibussowitsch               PetscCall(PetscHMapIGetWithDefault(cmap_i, bmap[cworkB[k]] + 1, 0, &tt));
1033a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1034a42ab0d6SHong Zhang             }
1035a42ab0d6SHong Zhang 
1036a42ab0d6SHong Zhang #else
1037a42ab0d6SHong Zhang             for (k = 0; k < nzA; k++) {
1038a42ab0d6SHong Zhang               if (cmap_i[cstart + cworkA[k]]) lens_i[j]++;
1039a42ab0d6SHong Zhang             }
1040a42ab0d6SHong Zhang             for (k = 0; k < nzB; k++) {
1041a42ab0d6SHong Zhang               if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++;
1042a42ab0d6SHong Zhang             }
1043a42ab0d6SHong Zhang #endif
1044a42ab0d6SHong Zhang           } else { /* allcolumns */
1045a42ab0d6SHong Zhang             lens_i[j] = nzA + nzB;
1046a42ab0d6SHong Zhang           }
104717df9f7cSHong Zhang         }
104817df9f7cSHong Zhang       }
104917df9f7cSHong Zhang     }
105017df9f7cSHong Zhang 
105117df9f7cSHong Zhang     /* Create row map: global row of C -> local row of submatrices */
105217df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
1053059f6b74SHong Zhang       if (!allrows[i]) {
1054059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE)
1055*eec179cfSJacob Faibussowitsch         PetscCall(PetscHMapICreateWithSize(nrow[i], rmap + i));
105617df9f7cSHong Zhang         irow_i = irow[i];
105717df9f7cSHong Zhang         jmax   = nrow[i];
105817df9f7cSHong Zhang         for (j = 0; j < jmax; j++) {
1059ec3c739cSHong Zhang           if (allrows[i]) {
1060*eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapISetWithMode(rmap[i], j + 1, j + 1, INSERT_VALUES));
1061ec3c739cSHong Zhang           } else {
1062*eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapISetWithMode(rmap[i], irow_i[j] + 1, j + 1, INSERT_VALUES));
106317df9f7cSHong Zhang           }
106417df9f7cSHong Zhang         }
106517df9f7cSHong Zhang #else
10669566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(c->Mbs, &rmap[i]));
106717df9f7cSHong Zhang         rmap_i = rmap[i];
106817df9f7cSHong Zhang         irow_i = irow[i];
106917df9f7cSHong Zhang         jmax   = nrow[i];
107017df9f7cSHong Zhang         for (j = 0; j < jmax; j++) {
1071ec3c739cSHong Zhang           if (allrows[i]) rmap_i[j] = j;
1072ec3c739cSHong Zhang           else rmap_i[irow_i[j]] = j;
107317df9f7cSHong Zhang         }
107417df9f7cSHong Zhang #endif
1075059f6b74SHong Zhang       } else rmap[i] = NULL;
1076059f6b74SHong Zhang     }
107717df9f7cSHong Zhang 
107817df9f7cSHong Zhang     /* Update lens from offproc data */
107917df9f7cSHong Zhang     {
108017df9f7cSHong Zhang       PetscInt *rbuf2_i, *rbuf3_i, *sbuf1_i;
108117df9f7cSHong Zhang 
10829566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Waitall(nrqs, r_waits3, MPI_STATUSES_IGNORE));
108317df9f7cSHong Zhang       for (tmp2 = 0; tmp2 < nrqs; tmp2++) {
108417df9f7cSHong Zhang         sbuf1_i = sbuf1[pa[tmp2]];
108517df9f7cSHong Zhang         jmax    = sbuf1_i[0];
108617df9f7cSHong Zhang         ct1     = 2 * jmax + 1;
108717df9f7cSHong Zhang         ct2     = 0;
108817df9f7cSHong Zhang         rbuf2_i = rbuf2[tmp2];
108917df9f7cSHong Zhang         rbuf3_i = rbuf3[tmp2];
109017df9f7cSHong Zhang         for (j = 1; j <= jmax; j++) {
109117df9f7cSHong Zhang           is_no  = sbuf1_i[2 * j - 1];
109217df9f7cSHong Zhang           max1   = sbuf1_i[2 * j];
109317df9f7cSHong Zhang           lens_i = lens[is_no];
109417df9f7cSHong Zhang           if (!allcolumns[is_no]) cmap_i = cmap[is_no];
109517df9f7cSHong Zhang           rmap_i = rmap[is_no];
109617df9f7cSHong Zhang           for (k = 0; k < max1; k++, ct1++) {
1097059f6b74SHong Zhang             if (allrows[is_no]) {
1098059f6b74SHong Zhang               row = sbuf1_i[ct1];
1099059f6b74SHong Zhang             } else {
110017df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1101*eec179cfSJacob Faibussowitsch               PetscCall(PetscHMapIGetWithDefault(rmap_i, sbuf1_i[ct1] + 1, 0, &row));
110217df9f7cSHong Zhang               row--;
110308401ef6SPierre Jolivet               PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table");
110417df9f7cSHong Zhang #else
110517df9f7cSHong Zhang               row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
110617df9f7cSHong Zhang #endif
1107059f6b74SHong Zhang             }
110817df9f7cSHong Zhang             max2 = rbuf2_i[ct1];
110917df9f7cSHong Zhang             for (l = 0; l < max2; l++, ct2++) {
111017df9f7cSHong Zhang               if (!allcolumns[is_no]) {
111117df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1112*eec179cfSJacob Faibussowitsch                 PetscCall(PetscHMapIGetWithDefault(cmap_i, rbuf3_i[ct2] + 1, 0, &tcol));
111317df9f7cSHong Zhang #else
111417df9f7cSHong Zhang                 tcol = cmap_i[rbuf3_i[ct2]];
111517df9f7cSHong Zhang #endif
111617df9f7cSHong Zhang                 if (tcol) lens_i[row]++;
111717df9f7cSHong Zhang               } else {         /* allcolumns */
111817df9f7cSHong Zhang                 lens_i[row]++; /* lens_i[row] += max2 ? */
111917df9f7cSHong Zhang               }
112017df9f7cSHong Zhang             }
112117df9f7cSHong Zhang           }
112217df9f7cSHong Zhang         }
112317df9f7cSHong Zhang       }
112417df9f7cSHong Zhang     }
11259566063dSJacob Faibussowitsch     PetscCall(PetscFree(r_waits3));
11269566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, s_waits3, MPI_STATUSES_IGNORE));
11279566063dSJacob Faibussowitsch     PetscCall(PetscFree(s_waits3));
112817df9f7cSHong Zhang 
112917df9f7cSHong Zhang     /* Create the submatrices */
113017df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
1131a42ab0d6SHong Zhang       PetscInt bs_tmp;
1132a42ab0d6SHong Zhang       if (ijonly) bs_tmp = 1;
1133a42ab0d6SHong Zhang       else bs_tmp = bs;
113417df9f7cSHong Zhang 
11359566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, submats + i));
11369566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(submats[i], nrow[i] * bs_tmp, ncol[i] * bs_tmp, PETSC_DETERMINE, PETSC_DETERMINE));
113717df9f7cSHong Zhang 
11389566063dSJacob Faibussowitsch       PetscCall(MatSetType(submats[i], ((PetscObject)A)->type_name));
11399566063dSJacob Faibussowitsch       PetscCall(MatSeqBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i]));
11409566063dSJacob Faibussowitsch       PetscCall(MatSeqSBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i])); /* this subroutine is used by SBAIJ routines */
114117df9f7cSHong Zhang 
11425c39f6d9SHong Zhang       /* create struct Mat_SubSppt and attached it to submat */
11439566063dSJacob Faibussowitsch       PetscCall(PetscNew(&smat_i));
114417df9f7cSHong Zhang       subc            = (Mat_SeqBAIJ *)submats[i]->data;
114517df9f7cSHong Zhang       subc->submatis1 = smat_i;
114617df9f7cSHong Zhang 
114717df9f7cSHong Zhang       smat_i->destroy          = submats[i]->ops->destroy;
1148f68bb481SHong Zhang       submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ;
114917df9f7cSHong Zhang       submats[i]->factortype   = C->factortype;
115017df9f7cSHong Zhang 
115117df9f7cSHong Zhang       smat_i->id          = i;
115217df9f7cSHong Zhang       smat_i->nrqs        = nrqs;
115317df9f7cSHong Zhang       smat_i->nrqr        = nrqr;
115417df9f7cSHong Zhang       smat_i->rbuf1       = rbuf1;
115517df9f7cSHong Zhang       smat_i->rbuf2       = rbuf2;
115617df9f7cSHong Zhang       smat_i->rbuf3       = rbuf3;
115717df9f7cSHong Zhang       smat_i->sbuf2       = sbuf2;
115817df9f7cSHong Zhang       smat_i->req_source2 = req_source2;
115917df9f7cSHong Zhang 
116017df9f7cSHong Zhang       smat_i->sbuf1 = sbuf1;
116117df9f7cSHong Zhang       smat_i->ptr   = ptr;
116217df9f7cSHong Zhang       smat_i->tmp   = tmp;
116317df9f7cSHong Zhang       smat_i->ctr   = ctr;
116417df9f7cSHong Zhang 
116517df9f7cSHong Zhang       smat_i->pa          = pa;
116617df9f7cSHong Zhang       smat_i->req_size    = req_size;
116717df9f7cSHong Zhang       smat_i->req_source1 = req_source1;
116817df9f7cSHong Zhang 
116917df9f7cSHong Zhang       smat_i->allcolumns = allcolumns[i];
1170ec3c739cSHong Zhang       smat_i->allrows    = allrows[i];
117117df9f7cSHong Zhang       smat_i->singleis   = PETSC_FALSE;
117217df9f7cSHong Zhang       smat_i->row2proc   = row2proc[i];
117317df9f7cSHong Zhang       smat_i->rmap       = rmap[i];
117417df9f7cSHong Zhang       smat_i->cmap       = cmap[i];
117517df9f7cSHong Zhang     }
117617df9f7cSHong Zhang 
11774698041cSHong Zhang     if (!ismax) { /* Create dummy submats[0] for reuse struct subc */
11789566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &submats[0]));
11799566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(submats[0], 0, 0, PETSC_DETERMINE, PETSC_DETERMINE));
11809566063dSJacob Faibussowitsch       PetscCall(MatSetType(submats[0], MATDUMMY));
11814698041cSHong Zhang 
11825c39f6d9SHong Zhang       /* create struct Mat_SubSppt and attached it to submat */
11834dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&smat_i));
11844698041cSHong Zhang       submats[0]->data = (void *)smat_i;
11854698041cSHong Zhang 
11864698041cSHong Zhang       smat_i->destroy          = submats[0]->ops->destroy;
1187f68bb481SHong Zhang       submats[0]->ops->destroy = MatDestroySubMatrix_Dummy;
11884698041cSHong Zhang       submats[0]->factortype   = C->factortype;
11894698041cSHong Zhang 
1190006cef31SHong Zhang       smat_i->id          = 0;
11914698041cSHong Zhang       smat_i->nrqs        = nrqs;
11924698041cSHong Zhang       smat_i->nrqr        = nrqr;
11934698041cSHong Zhang       smat_i->rbuf1       = rbuf1;
11944698041cSHong Zhang       smat_i->rbuf2       = rbuf2;
11954698041cSHong Zhang       smat_i->rbuf3       = rbuf3;
11964698041cSHong Zhang       smat_i->sbuf2       = sbuf2;
11974698041cSHong Zhang       smat_i->req_source2 = req_source2;
11984698041cSHong Zhang 
11994698041cSHong Zhang       smat_i->sbuf1 = sbuf1;
12004698041cSHong Zhang       smat_i->ptr   = ptr;
12014698041cSHong Zhang       smat_i->tmp   = tmp;
12024698041cSHong Zhang       smat_i->ctr   = ctr;
12034698041cSHong Zhang 
12044698041cSHong Zhang       smat_i->pa          = pa;
12054698041cSHong Zhang       smat_i->req_size    = req_size;
12064698041cSHong Zhang       smat_i->req_source1 = req_source1;
12074698041cSHong Zhang 
1208c67fe082SHong Zhang       smat_i->allcolumns = PETSC_FALSE;
12094698041cSHong Zhang       smat_i->singleis   = PETSC_FALSE;
12104698041cSHong Zhang       smat_i->row2proc   = NULL;
12114698041cSHong Zhang       smat_i->rmap       = NULL;
12124698041cSHong Zhang       smat_i->cmap       = NULL;
12134698041cSHong Zhang     }
12144698041cSHong Zhang 
12159566063dSJacob Faibussowitsch     if (ismax) PetscCall(PetscFree(lens[0]));
12169566063dSJacob Faibussowitsch     PetscCall(PetscFree(lens));
1217175c2bc5SJunchao Zhang     if (sbuf_aj) {
12189566063dSJacob Faibussowitsch       PetscCall(PetscFree(sbuf_aj[0]));
12199566063dSJacob Faibussowitsch       PetscCall(PetscFree(sbuf_aj));
1220175c2bc5SJunchao Zhang     }
122117df9f7cSHong Zhang 
122217df9f7cSHong Zhang   } /* endof scall == MAT_INITIAL_MATRIX */
122317df9f7cSHong Zhang 
122417df9f7cSHong Zhang   /* Post recv matrix values */
1225bbc89d27SHong Zhang   if (!ijonly) {
12269566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag4));
12279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &rbuf4));
12289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &r_waits4));
122917df9f7cSHong Zhang     for (i = 0; i < nrqs; ++i) {
12309566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(rbuf2[i][0] * bs2, &rbuf4[i]));
12319566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(rbuf4[i], rbuf2[i][0] * bs2, MPIU_SCALAR, req_source2[i], tag4, comm, r_waits4 + i));
123217df9f7cSHong Zhang     }
123317df9f7cSHong Zhang 
123417df9f7cSHong Zhang     /* Allocate sending buffers for a->a, and send them off */
12359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &sbuf_aa));
123617df9f7cSHong Zhang     for (i = 0, j = 0; i < nrqr; i++) j += req_size[i];
12379e686edcSHong Zhang 
12389566063dSJacob Faibussowitsch     if (nrqr) PetscCall(PetscMalloc1(j * bs2, &sbuf_aa[0]));
1239a42ab0d6SHong Zhang     for (i = 1; i < nrqr; i++) sbuf_aa[i] = sbuf_aa[i - 1] + req_size[i - 1] * bs2;
124017df9f7cSHong Zhang 
12419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &s_waits4));
124217df9f7cSHong Zhang 
124317df9f7cSHong Zhang     for (i = 0; i < nrqr; i++) {
124417df9f7cSHong Zhang       rbuf1_i   = rbuf1[i];
124517df9f7cSHong Zhang       sbuf_aa_i = sbuf_aa[i];
124617df9f7cSHong Zhang       ct1       = 2 * rbuf1_i[0] + 1;
124717df9f7cSHong Zhang       ct2       = 0;
124817df9f7cSHong Zhang       for (j = 1, max1 = rbuf1_i[0]; j <= max1; j++) {
124917df9f7cSHong Zhang         kmax = rbuf1_i[2 * j];
125017df9f7cSHong Zhang         for (k = 0; k < kmax; k++, ct1++) {
125117df9f7cSHong Zhang           row    = rbuf1_i[ct1] - rstart;
1252a42ab0d6SHong Zhang           nzA    = a_i[row + 1] - a_i[row];
1253a42ab0d6SHong Zhang           nzB    = b_i[row + 1] - b_i[row];
125417df9f7cSHong Zhang           ncols  = nzA + nzB;
125517df9f7cSHong Zhang           cworkB = b_j + b_i[row];
1256a42ab0d6SHong Zhang           vworkA = a_a + a_i[row] * bs2;
1257a42ab0d6SHong Zhang           vworkB = b_a + b_i[row] * bs2;
125817df9f7cSHong Zhang 
125917df9f7cSHong Zhang           /* load the column values for this row into vals*/
1260a42ab0d6SHong Zhang           vals = sbuf_aa_i + ct2 * bs2;
1261a42ab0d6SHong Zhang           for (l = 0; l < nzB; l++) {
1262a42ab0d6SHong Zhang             if ((bmap[cworkB[l]]) < cstart) {
12639566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(vals + l * bs2, vworkB + l * bs2, bs2));
1264a42ab0d6SHong Zhang             } else break;
1265a42ab0d6SHong Zhang           }
1266a42ab0d6SHong Zhang           imark = l;
126748a46eb9SPierre Jolivet           for (l = 0; l < nzA; l++) PetscCall(PetscArraycpy(vals + (imark + l) * bs2, vworkA + l * bs2, bs2));
126848a46eb9SPierre Jolivet           for (l = imark; l < nzB; l++) PetscCall(PetscArraycpy(vals + (nzA + l) * bs2, vworkB + l * bs2, bs2));
126917df9f7cSHong Zhang 
127017df9f7cSHong Zhang           ct2 += ncols;
127117df9f7cSHong Zhang         }
127217df9f7cSHong Zhang       }
12739566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sbuf_aa_i, req_size[i] * bs2, MPIU_SCALAR, req_source1[i], tag4, comm, s_waits4 + i));
127417df9f7cSHong Zhang     }
1275bbc89d27SHong Zhang   }
127617df9f7cSHong Zhang 
127717df9f7cSHong Zhang   /* Assemble the matrices */
127817df9f7cSHong Zhang   /* First assemble the local rows */
127917df9f7cSHong Zhang   for (i = 0; i < ismax; i++) {
128017df9f7cSHong Zhang     row2proc_i = row2proc[i];
128117df9f7cSHong Zhang     subc       = (Mat_SeqBAIJ *)submats[i]->data;
128217df9f7cSHong Zhang     imat_ilen  = subc->ilen;
128317df9f7cSHong Zhang     imat_j     = subc->j;
128417df9f7cSHong Zhang     imat_i     = subc->i;
128517df9f7cSHong Zhang     imat_a     = subc->a;
128617df9f7cSHong Zhang 
128717df9f7cSHong Zhang     if (!allcolumns[i]) cmap_i = cmap[i];
128817df9f7cSHong Zhang     rmap_i = rmap[i];
128917df9f7cSHong Zhang     irow_i = irow[i];
129017df9f7cSHong Zhang     jmax   = nrow[i];
129117df9f7cSHong Zhang     for (j = 0; j < jmax; j++) {
1292ec3c739cSHong Zhang       if (allrows[i]) row = j;
1293ec3c739cSHong Zhang       else row = irow_i[j];
129417df9f7cSHong Zhang       proc = row2proc_i[j];
1295a42ab0d6SHong Zhang 
129617df9f7cSHong Zhang       if (proc == rank) {
1297a42ab0d6SHong Zhang         row    = row - rstart;
1298a42ab0d6SHong Zhang         nzA    = a_i[row + 1] - a_i[row];
1299a42ab0d6SHong Zhang         nzB    = b_i[row + 1] - b_i[row];
1300a42ab0d6SHong Zhang         cworkA = a_j + a_i[row];
1301a42ab0d6SHong Zhang         cworkB = b_j + b_i[row];
1302a42ab0d6SHong Zhang         if (!ijonly) {
1303a42ab0d6SHong Zhang           vworkA = a_a + a_i[row] * bs2;
1304a42ab0d6SHong Zhang           vworkB = b_a + b_i[row] * bs2;
1305a42ab0d6SHong Zhang         }
1306059f6b74SHong Zhang 
1307059f6b74SHong Zhang         if (allrows[i]) {
1308059f6b74SHong Zhang           row = row + rstart;
1309059f6b74SHong Zhang         } else {
1310a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1311*eec179cfSJacob Faibussowitsch           PetscCall(PetscHMapIGetWithDefault(rmap_i, row + rstart + 1, 0, &row));
1312a42ab0d6SHong Zhang           row--;
1313121971b2SHong Zhang 
131408401ef6SPierre Jolivet           PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table");
1315a42ab0d6SHong Zhang #else
1316a42ab0d6SHong Zhang           row = rmap_i[row + rstart];
1317a42ab0d6SHong Zhang #endif
1318059f6b74SHong Zhang         }
1319a42ab0d6SHong Zhang         mat_i = imat_i[row];
1320a42ab0d6SHong Zhang         if (!ijonly) mat_a = imat_a + mat_i * bs2;
1321a42ab0d6SHong Zhang         mat_j = imat_j + mat_i;
132280d31651SHong Zhang         ilen  = imat_ilen[row];
1323a42ab0d6SHong Zhang 
1324a42ab0d6SHong Zhang         /* load the column indices for this row into cols*/
1325a42ab0d6SHong Zhang         if (!allcolumns[i]) {
1326a42ab0d6SHong Zhang           for (l = 0; l < nzB; l++) {
1327a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1328a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1329*eec179cfSJacob Faibussowitsch               PetscCall(PetscHMapIGetWithDefault(cmap_i, ctmp + 1, 0, &tcol));
1330a42ab0d6SHong Zhang               if (tcol) {
1331a42ab0d6SHong Zhang #else
1332a42ab0d6SHong Zhang               if ((tcol = cmap_i[ctmp])) {
1333a42ab0d6SHong Zhang #endif
1334a42ab0d6SHong Zhang                 *mat_j++ = tcol - 1;
13359566063dSJacob Faibussowitsch                 PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2));
1336a42ab0d6SHong Zhang                 mat_a += bs2;
133780d31651SHong Zhang                 ilen++;
1338a42ab0d6SHong Zhang               }
1339a42ab0d6SHong Zhang             } else break;
1340a42ab0d6SHong Zhang           }
1341a42ab0d6SHong Zhang           imark = l;
1342a42ab0d6SHong Zhang           for (l = 0; l < nzA; l++) {
1343a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1344*eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapIGetWithDefault(cmap_i, cstart + cworkA[l] + 1, 0, &tcol));
1345a42ab0d6SHong Zhang             if (tcol) {
1346a42ab0d6SHong Zhang #else
1347a42ab0d6SHong Zhang             if ((tcol = cmap_i[cstart + cworkA[l]])) {
1348a42ab0d6SHong Zhang #endif
1349a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1350a42ab0d6SHong Zhang               if (!ijonly) {
13519566063dSJacob Faibussowitsch                 PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2));
1352a42ab0d6SHong Zhang                 mat_a += bs2;
1353a42ab0d6SHong Zhang               }
135480d31651SHong Zhang               ilen++;
1355a42ab0d6SHong Zhang             }
1356a42ab0d6SHong Zhang           }
1357a42ab0d6SHong Zhang           for (l = imark; l < nzB; l++) {
1358a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1359*eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapIGetWithDefault(cmap_i, bmap[cworkB[l]] + 1, 0, &tcol));
1360a42ab0d6SHong Zhang             if (tcol) {
1361a42ab0d6SHong Zhang #else
1362a42ab0d6SHong Zhang             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1363a42ab0d6SHong Zhang #endif
1364a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1365a42ab0d6SHong Zhang               if (!ijonly) {
13669566063dSJacob Faibussowitsch                 PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2));
1367a42ab0d6SHong Zhang                 mat_a += bs2;
1368a42ab0d6SHong Zhang               }
136980d31651SHong Zhang               ilen++;
1370a42ab0d6SHong Zhang             }
1371a42ab0d6SHong Zhang           }
1372a42ab0d6SHong Zhang         } else { /* allcolumns */
1373a42ab0d6SHong Zhang           for (l = 0; l < nzB; l++) {
1374a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1375a42ab0d6SHong Zhang               *mat_j++ = ctmp;
13769566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2));
1377a42ab0d6SHong Zhang               mat_a += bs2;
137880d31651SHong Zhang               ilen++;
1379a42ab0d6SHong Zhang             } else break;
1380a42ab0d6SHong Zhang           }
1381a42ab0d6SHong Zhang           imark = l;
1382a42ab0d6SHong Zhang           for (l = 0; l < nzA; l++) {
1383a42ab0d6SHong Zhang             *mat_j++ = cstart + cworkA[l];
1384a42ab0d6SHong Zhang             if (!ijonly) {
13859566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2));
1386a42ab0d6SHong Zhang               mat_a += bs2;
1387a42ab0d6SHong Zhang             }
138880d31651SHong Zhang             ilen++;
1389a42ab0d6SHong Zhang           }
1390a42ab0d6SHong Zhang           for (l = imark; l < nzB; l++) {
1391a42ab0d6SHong Zhang             *mat_j++ = bmap[cworkB[l]];
1392a42ab0d6SHong Zhang             if (!ijonly) {
13939566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2));
1394a42ab0d6SHong Zhang               mat_a += bs2;
1395a42ab0d6SHong Zhang             }
139680d31651SHong Zhang             ilen++;
1397a42ab0d6SHong Zhang           }
1398a42ab0d6SHong Zhang         }
139980d31651SHong Zhang         imat_ilen[row] = ilen;
140017df9f7cSHong Zhang       }
140117df9f7cSHong Zhang     }
140217df9f7cSHong Zhang   }
140317df9f7cSHong Zhang 
140417df9f7cSHong Zhang   /* Now assemble the off proc rows */
140548a46eb9SPierre Jolivet   if (!ijonly) PetscCallMPI(MPI_Waitall(nrqs, r_waits4, MPI_STATUSES_IGNORE));
140617df9f7cSHong Zhang   for (tmp2 = 0; tmp2 < nrqs; tmp2++) {
140717df9f7cSHong Zhang     sbuf1_i = sbuf1[pa[tmp2]];
140817df9f7cSHong Zhang     jmax    = sbuf1_i[0];
140917df9f7cSHong Zhang     ct1     = 2 * jmax + 1;
141017df9f7cSHong Zhang     ct2     = 0;
141117df9f7cSHong Zhang     rbuf2_i = rbuf2[tmp2];
141217df9f7cSHong Zhang     rbuf3_i = rbuf3[tmp2];
14135dd0c08cSHong Zhang     if (!ijonly) rbuf4_i = rbuf4[tmp2];
141417df9f7cSHong Zhang     for (j = 1; j <= jmax; j++) {
141517df9f7cSHong Zhang       is_no  = sbuf1_i[2 * j - 1];
141617df9f7cSHong Zhang       rmap_i = rmap[is_no];
141717df9f7cSHong Zhang       if (!allcolumns[is_no]) cmap_i = cmap[is_no];
141817df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ *)submats[is_no]->data;
141917df9f7cSHong Zhang       imat_ilen = subc->ilen;
142017df9f7cSHong Zhang       imat_j    = subc->j;
142117df9f7cSHong Zhang       imat_i    = subc->i;
14225dd0c08cSHong Zhang       if (!ijonly) imat_a = subc->a;
142317df9f7cSHong Zhang       max1 = sbuf1_i[2 * j];
14245dd0c08cSHong Zhang       for (k = 0; k < max1; k++, ct1++) { /* for each recved block row */
142517df9f7cSHong Zhang         row = sbuf1_i[ct1];
1426059f6b74SHong Zhang 
1427059f6b74SHong Zhang         if (allrows[is_no]) {
1428059f6b74SHong Zhang           row = sbuf1_i[ct1];
1429059f6b74SHong Zhang         } else {
143017df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1431*eec179cfSJacob Faibussowitsch           PetscCall(PetscHMapIGetWithDefault(rmap_i, row + 1, 0, &row));
143217df9f7cSHong Zhang           row--;
143308401ef6SPierre Jolivet           PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table");
143417df9f7cSHong Zhang #else
143517df9f7cSHong Zhang           row = rmap_i[row];
143617df9f7cSHong Zhang #endif
1437059f6b74SHong Zhang         }
143817df9f7cSHong Zhang         ilen  = imat_ilen[row];
143917df9f7cSHong Zhang         mat_i = imat_i[row];
14405dd0c08cSHong Zhang         if (!ijonly) mat_a = imat_a + mat_i * bs2;
144117df9f7cSHong Zhang         mat_j = imat_j + mat_i;
144217df9f7cSHong Zhang         max2  = rbuf2_i[ct1];
144317df9f7cSHong Zhang         if (!allcolumns[is_no]) {
144417df9f7cSHong Zhang           for (l = 0; l < max2; l++, ct2++) {
144517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1446*eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapIGetWithDefault(cmap_i, rbuf3_i[ct2] + 1, 0, &tcol));
144717df9f7cSHong Zhang #else
144817df9f7cSHong Zhang             tcol = cmap_i[rbuf3_i[ct2]];
144917df9f7cSHong Zhang #endif
145017df9f7cSHong Zhang             if (tcol) {
145117df9f7cSHong Zhang               *mat_j++ = tcol - 1;
14525dd0c08cSHong Zhang               if (!ijonly) {
14539566063dSJacob Faibussowitsch                 PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2));
14545dd0c08cSHong Zhang                 mat_a += bs2;
14555dd0c08cSHong Zhang               }
145617df9f7cSHong Zhang               ilen++;
145717df9f7cSHong Zhang             }
145817df9f7cSHong Zhang           }
145917df9f7cSHong Zhang         } else { /* allcolumns */
146017df9f7cSHong Zhang           for (l = 0; l < max2; l++, ct2++) {
146117df9f7cSHong Zhang             *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */
14625dd0c08cSHong Zhang             if (!ijonly) {
14639566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2));
14645dd0c08cSHong Zhang               mat_a += bs2;
14655dd0c08cSHong Zhang             }
146617df9f7cSHong Zhang             ilen++;
146717df9f7cSHong Zhang           }
146817df9f7cSHong Zhang         }
146917df9f7cSHong Zhang         imat_ilen[row] = ilen;
147017df9f7cSHong Zhang       }
147117df9f7cSHong Zhang     }
147217df9f7cSHong Zhang   }
147317df9f7cSHong Zhang 
147480d31651SHong Zhang   if (!iscsorted) { /* sort column indices of the rows */
147580d31651SHong Zhang     MatScalar *work;
147680d31651SHong Zhang 
14779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs2, &work));
147817df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
147917df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ *)submats[i]->data;
148080d31651SHong Zhang       imat_ilen = subc->ilen;
148117df9f7cSHong Zhang       imat_j    = subc->j;
148217df9f7cSHong Zhang       imat_i    = subc->i;
14834b6ceb0dSHong Zhang       if (!ijonly) imat_a = subc->a;
148417df9f7cSHong Zhang       if (allcolumns[i]) continue;
14854b6ceb0dSHong Zhang 
148617df9f7cSHong Zhang       jmax = nrow[i];
148717df9f7cSHong Zhang       for (j = 0; j < jmax; j++) {
148880d31651SHong Zhang         mat_i = imat_i[j];
148980d31651SHong Zhang         mat_j = imat_j + mat_i;
14904b6ceb0dSHong Zhang         ilen  = imat_ilen[j];
149180d31651SHong Zhang         if (ijonly) {
14929566063dSJacob Faibussowitsch           PetscCall(PetscSortInt(ilen, mat_j));
149380d31651SHong Zhang         } else {
14944b6ceb0dSHong Zhang           mat_a = imat_a + mat_i * bs2;
14959566063dSJacob Faibussowitsch           PetscCall(PetscSortIntWithDataArray(ilen, mat_j, mat_a, bs2 * sizeof(MatScalar), work));
149617df9f7cSHong Zhang         }
149717df9f7cSHong Zhang       }
149817df9f7cSHong Zhang     }
14999566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
150080d31651SHong Zhang   }
150117df9f7cSHong Zhang 
1502bbc89d27SHong Zhang   if (!ijonly) {
15039566063dSJacob Faibussowitsch     PetscCall(PetscFree(r_waits4));
15049566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, s_waits4, MPI_STATUSES_IGNORE));
15059566063dSJacob Faibussowitsch     PetscCall(PetscFree(s_waits4));
1506bbc89d27SHong Zhang   }
150717df9f7cSHong Zhang 
150817df9f7cSHong Zhang   /* Restore the indices */
150917df9f7cSHong Zhang   for (i = 0; i < ismax; i++) {
151048a46eb9SPierre Jolivet     if (!allrows[i]) PetscCall(ISRestoreIndices(isrow[i], irow + i));
151148a46eb9SPierre Jolivet     if (!allcolumns[i]) PetscCall(ISRestoreIndices(iscol[i], icol + i));
151217df9f7cSHong Zhang   }
151317df9f7cSHong Zhang 
151417df9f7cSHong Zhang   for (i = 0; i < ismax; i++) {
15159566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY));
15169566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY));
151717df9f7cSHong Zhang   }
151817df9f7cSHong Zhang 
15199566063dSJacob Faibussowitsch   PetscCall(PetscFree5(*(PetscInt ***)&irow, *(PetscInt ***)&icol, nrow, ncol, issorted));
15209566063dSJacob Faibussowitsch   PetscCall(PetscFree5(row2proc, cmap, rmap, allcolumns, allrows));
152117df9f7cSHong Zhang 
1522bbc89d27SHong Zhang   if (!ijonly) {
1523175c2bc5SJunchao Zhang     if (sbuf_aa) {
15249566063dSJacob Faibussowitsch       PetscCall(PetscFree(sbuf_aa[0]));
15259566063dSJacob Faibussowitsch       PetscCall(PetscFree(sbuf_aa));
1526175c2bc5SJunchao Zhang     }
152717df9f7cSHong Zhang 
152848a46eb9SPierre Jolivet     for (i = 0; i < nrqs; ++i) PetscCall(PetscFree(rbuf4[i]));
15299566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf4));
1530bbc89d27SHong Zhang   }
15317a868f3eSHong Zhang   c->ijonly = PETSC_FALSE; /* set back to the default */
15323a40ed3dSBarry Smith   PetscFunctionReturn(0);
1533d5b485abSSatish Balay }
1534