xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
1d5b485abSSatish Balay /*
2d5b485abSSatish Balay    Routines to compute overlapping regions of a parallel MPI matrix
3d5b485abSSatish Balay   and to find submatrices that were shared across processors.
4d5b485abSSatish Balay */
5c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h>
6c6db04a5SJed Brown #include <petscbt.h>
7d5b485abSSatish Balay 
8b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat, PetscInt, char **, PetscInt *, PetscInt **);
9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat, PetscInt, PetscInt **, PetscInt **, PetscInt *);
1009573ac7SBarry Smith extern PetscErrorCode MatGetRow_MPIBAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
1109573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
12d5b485abSSatish Balay 
13d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C, PetscInt imax, IS is[], PetscInt ov)
14d71ae5a4SJacob Faibussowitsch {
15bd38d1d0SPierre Jolivet   PetscInt        i, N = C->cmap->N, bs = C->rmap->bs, n;
16bd38d1d0SPierre Jolivet   const PetscInt *idx;
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));
25bd38d1d0SPierre Jolivet   for (i = 0; i < imax; i++) {
26bd38d1d0SPierre Jolivet     PetscCall(ISDestroy(&is[i]));
27bd38d1d0SPierre Jolivet     PetscCall(ISGetLocalSize(is_new[i], &n));
28bd38d1d0SPierre Jolivet     PetscCall(ISGetIndices(is_new[i], &idx));
29bd38d1d0SPierre Jolivet     PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)is_new[i]), bs, n, idx, PETSC_COPY_VALUES, &is[i]));
30bd38d1d0SPierre Jolivet     PetscCall(ISDestroy(&is_new[i]));
31bd38d1d0SPierre Jolivet   }
329566063dSJacob Faibussowitsch   PetscCall(PetscFree(is_new));
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34d5b485abSSatish Balay }
35d5b485abSSatish Balay 
36d5b485abSSatish Balay /*
37d5b485abSSatish Balay   Sample message format:
38d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
390e9b0e7eSHong Zhang   to index sets is[1], is[5]
40d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
41d5b485abSSatish Balay   -----------
42d5b485abSSatish Balay   mesg [1] = 1 => is[1]
43d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
44d5b485abSSatish Balay   -----------
45d5b485abSSatish Balay   mesg [5] = 5  => is[5]
46d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
47d5b485abSSatish Balay   -----------
48d5b485abSSatish Balay   mesg [7]
490e9b0e7eSHong Zhang   mesg [n]  data(is[1])
50d5b485abSSatish Balay   -----------
51d5b485abSSatish Balay   mesg[n+1]
52d5b485abSSatish Balay   mesg[m]  data(is[5])
53d5b485abSSatish Balay   -----------
54d5b485abSSatish Balay 
55d5b485abSSatish Balay   Notes:
56d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
572d4ee042Sprj-   nrqr - no of requests received (which have to be or which have been processed)
58d5b485abSSatish Balay */
59d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C, PetscInt imax, IS is[])
60d71ae5a4SJacob Faibussowitsch {
61df36731bSBarry Smith   Mat_MPIBAIJ     *c = (Mat_MPIBAIJ *)C->data;
625d0c19d7SBarry Smith   const PetscInt **idx, *idx_i;
6324c049a4SHong Zhang   PetscInt        *n, *w3, *w4, **data, len;
64b24ad042SBarry Smith   PetscMPIInt      size, rank, tag1, tag2, *w2, *w1, nrqr;
65131c27b5Sprj-   PetscInt         Mbs, i, j, k, **rbuf, row, nrqs, msz, **outdat, **ptr;
668f9f447aSHong Zhang   PetscInt        *ctr, *pa, *tmp, *isz, *isz1, **xdata, **rbuf2, *d_p;
67131c27b5Sprj-   PetscMPIInt     *onodes1, *olengths1, *onodes2, *olengths2, proc = -1;
68f1af5d2fSBarry Smith   PetscBT         *table;
69749130a8SStefano Zampini   MPI_Comm         comm, *iscomms;
70d5b485abSSatish Balay   MPI_Request     *s_waits1, *r_waits1, *s_waits2, *r_waits2;
718f9f447aSHong Zhang   char            *t_p;
72d5b485abSSatish Balay 
733a40ed3dSBarry Smith   PetscFunctionBegin;
749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
75d5b485abSSatish Balay   size = c->size;
76d5b485abSSatish Balay   rank = c->rank;
77df36731bSBarry Smith   Mbs  = c->Mbs;
78d5b485abSSatish Balay 
799566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag1));
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2));
81c7dd2924SSatish Balay 
82864bd19bSPierre Jolivet   PetscCall(PetscMalloc2(imax, (PetscInt ***)&idx, imax, &n));
8324c049a4SHong Zhang 
84d5b485abSSatish Balay   for (i = 0; i < imax; i++) {
859566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[i], &idx[i]));
869566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(is[i], &n[i]));
87d5b485abSSatish Balay   }
88d5b485abSSatish Balay 
89d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
90d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
919566063dSJacob Faibussowitsch   PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4));
92d5b485abSSatish Balay   for (i = 0; i < imax; i++) {
939566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(w4, size)); /* initialise work vector*/
94d5b485abSSatish Balay     idx_i = idx[i];
95d5b485abSSatish Balay     len   = n[i];
96d5b485abSSatish Balay     for (j = 0; j < len; j++) {
97d5b485abSSatish Balay       row = idx_i[j];
9808401ef6SPierre Jolivet       PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index set cannot have negative entries");
999566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc));
100d5b485abSSatish Balay       w4[proc]++;
101d5b485abSSatish Balay     }
102d5b485abSSatish Balay     for (j = 0; j < size; j++) {
1039371c9d4SSatish Balay       if (w4[j]) {
1049371c9d4SSatish Balay         w1[j] += w4[j];
1059371c9d4SSatish Balay         w3[j]++;
1069371c9d4SSatish Balay       }
107d5b485abSSatish Balay     }
108d5b485abSSatish Balay   }
109d5b485abSSatish Balay 
110d5b485abSSatish Balay   nrqs     = 0; /* no of outgoing messages */
111d5b485abSSatish Balay   msz      = 0; /* total mesg length (for all proc */
1120e9b0e7eSHong Zhang   w1[rank] = 0; /* no mesg sent to itself */
113d5b485abSSatish Balay   w3[rank] = 0;
114d5b485abSSatish Balay   for (i = 0; i < size; i++) {
1159371c9d4SSatish Balay     if (w1[i]) {
1169371c9d4SSatish Balay       w2[i] = 1;
1179371c9d4SSatish Balay       nrqs++;
1189371c9d4SSatish Balay     } /* there exists a message to proc i */
119d5b485abSSatish Balay   }
120d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
1219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqs, &pa));
122d5b485abSSatish Balay   for (i = 0, j = 0; i < size; i++) {
1239371c9d4SSatish Balay     if (w1[i]) {
1249371c9d4SSatish Balay       pa[j] = i;
1259371c9d4SSatish Balay       j++;
1269371c9d4SSatish Balay     }
127d5b485abSSatish Balay   }
128d5b485abSSatish Balay 
129d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
130d5b485abSSatish Balay   for (i = 0; i < nrqs; i++) {
131d5b485abSSatish Balay     j = pa[i];
132d5b485abSSatish Balay     w1[j] += w2[j] + 2 * w3[j];
133d5b485abSSatish Balay     msz += w1[j];
134d5b485abSSatish Balay   }
135d5b485abSSatish Balay 
136c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
1379566063dSJacob Faibussowitsch   PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr));
1389566063dSJacob Faibussowitsch   PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1));
139d5b485abSSatish Balay 
140c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
1419566063dSJacob Faibussowitsch   PetscCall(PetscPostIrecvInt(comm, tag1, nrqr, onodes1, olengths1, &rbuf, &r_waits1));
142d5b485abSSatish Balay 
143d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
1449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(size, &outdat, size, &ptr, msz, &tmp, size, &ctr));
1459566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(outdat, size));
1469566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(ptr, size));
147d5b485abSSatish Balay   {
148b24ad042SBarry Smith     PetscInt *iptr = tmp, ict = 0;
149d5b485abSSatish Balay     for (i = 0; i < nrqs; i++) {
150d5b485abSSatish Balay       j = pa[i];
151d5b485abSSatish Balay       iptr += ict;
152d5b485abSSatish Balay       outdat[j] = iptr;
153d5b485abSSatish Balay       ict       = w1[j];
154d5b485abSSatish Balay     }
155d5b485abSSatish Balay   }
156d5b485abSSatish Balay 
157d5b485abSSatish Balay   /* Form the outgoing messages */
158d5b485abSSatish Balay   /*plug in the headers*/
159d5b485abSSatish Balay   for (i = 0; i < nrqs; i++) {
160d5b485abSSatish Balay     j            = pa[i];
161d5b485abSSatish Balay     outdat[j][0] = 0;
1629566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(outdat[j] + 1, 2 * w3[j]));
163d5b485abSSatish Balay     ptr[j] = outdat[j] + 2 * w3[j] + 1;
164d5b485abSSatish Balay   }
165d5b485abSSatish Balay 
166d5b485abSSatish Balay   /* Memory for doing local proc's work*/
167d5b485abSSatish Balay   {
1689566063dSJacob Faibussowitsch     PetscCall(PetscCalloc5(imax, &table, imax, &data, imax, &isz, Mbs * imax, &d_p, (Mbs / PETSC_BITS_PER_BYTE + 1) * imax, &t_p));
169d5b485abSSatish Balay 
170bbd702dbSSatish Balay     for (i = 0; i < imax; i++) {
171b6410449SSatish Balay       table[i] = t_p + (Mbs / PETSC_BITS_PER_BYTE + 1) * i;
172bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
173d5b485abSSatish Balay     }
174d5b485abSSatish Balay   }
175d5b485abSSatish Balay 
176d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
177d5b485abSSatish Balay   {
178b24ad042SBarry Smith     PetscInt n_i, *data_i, isz_i, *outdat_j, ctr_j;
179f1af5d2fSBarry Smith     PetscBT  table_i;
180d5b485abSSatish Balay 
181d5b485abSSatish Balay     for (i = 0; i < imax; i++) {
1829566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctr, size));
183d5b485abSSatish Balay       n_i     = n[i];
184d5b485abSSatish Balay       table_i = table[i];
185d5b485abSSatish Balay       idx_i   = idx[i];
186d5b485abSSatish Balay       data_i  = data[i];
187d5b485abSSatish Balay       isz_i   = isz[i];
188d5b485abSSatish Balay       for (j = 0; j < n_i; j++) { /* parse the indices of each IS */
189d5b485abSSatish Balay         row = idx_i[j];
1909566063dSJacob Faibussowitsch         PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc));
191d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
192d5b485abSSatish Balay           ctr[proc]++;
193d5b485abSSatish Balay           *ptr[proc] = row;
194d5b485abSSatish Balay           ptr[proc]++;
195d6b45a43SBarry Smith         } else { /* Update the local table */
19626fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row;
197d5b485abSSatish Balay         }
198d5b485abSSatish Balay       }
199d5b485abSSatish Balay       /* Update the headers for the current IS */
200d5b485abSSatish Balay       for (j = 0; j < size; j++) { /* Can Optimise this loop by using pa[] */
201d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
202d5b485abSSatish Balay           outdat_j            = outdat[j];
203d5b485abSSatish Balay           k                   = ++outdat_j[0];
204d5b485abSSatish Balay           outdat_j[2 * k]     = ctr_j;
205d5b485abSSatish Balay           outdat_j[2 * k - 1] = i;
206d5b485abSSatish Balay         }
207d5b485abSSatish Balay       }
208d5b485abSSatish Balay       isz[i] = isz_i;
209d5b485abSSatish Balay     }
210d5b485abSSatish Balay   }
211d5b485abSSatish Balay 
212d5b485abSSatish Balay   /*  Now  post the sends */
2139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqs, &s_waits1));
214d5b485abSSatish Balay   for (i = 0; i < nrqs; ++i) {
215d5b485abSSatish Balay     j = pa[i];
2169566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Isend(outdat[j], w1[j], MPIU_INT, j, tag1, comm, s_waits1 + i));
217d5b485abSSatish Balay   }
218d5b485abSSatish Balay 
219d5b485abSSatish Balay   /* No longer need the original indices*/
22048a46eb9SPierre Jolivet   for (i = 0; i < imax; ++i) PetscCall(ISRestoreIndices(is[i], idx + i));
2219566063dSJacob Faibussowitsch   PetscCall(PetscFree2(*(PetscInt ***)&idx, n));
222d5b485abSSatish Balay 
2239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(imax, &iscomms));
224d5b485abSSatish Balay   for (i = 0; i < imax; ++i) {
2259566063dSJacob Faibussowitsch     PetscCall(PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]), &iscomms[i], NULL));
2269566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is[i]));
227d5b485abSSatish Balay   }
228d5b485abSSatish Balay 
229d5b485abSSatish Balay   /* Do Local work*/
2309566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap_MPIBAIJ_Local(C, imax, table, isz, data));
231d5b485abSSatish Balay 
232d5b485abSSatish Balay   /* Receive messages*/
2339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE));
2349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE));
235d5b485abSSatish Balay 
236d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
2379566063dSJacob Faibussowitsch   PetscCall(PetscFree4(outdat, ptr, tmp, ctr));
2389566063dSJacob Faibussowitsch   PetscCall(PetscFree4(w1, w2, w3, w4));
239d5b485abSSatish Balay 
2409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqr, &xdata));
2419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqr, &isz1));
2429566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap_MPIBAIJ_Receive(C, nrqr, rbuf, xdata, isz1));
243175c2bc5SJunchao Zhang   if (rbuf) {
2449566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf[0]));
2459566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf));
246175c2bc5SJunchao Zhang   }
247d5b485abSSatish Balay 
248d5b485abSSatish Balay   /* Send the data back*/
249d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
250d5b485abSSatish Balay   {
251b24ad042SBarry Smith     PetscMPIInt *rw1;
252d5b485abSSatish Balay 
2539566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(size, &rw1));
254c7dd2924SSatish Balay 
255d5b485abSSatish Balay     for (i = 0; i < nrqr; ++i) {
256175c2bc5SJunchao Zhang       proc      = onodes1[i];
257d5b485abSSatish Balay       rw1[proc] = isz1[i];
258d5b485abSSatish Balay     }
259d5b485abSSatish Balay 
260c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
2619566063dSJacob Faibussowitsch     PetscCall(PetscGatherMessageLengths(comm, nrqr, nrqs, rw1, &onodes2, &olengths2));
2629566063dSJacob Faibussowitsch     PetscCall(PetscFree(rw1));
263c7dd2924SSatish Balay   }
264c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
2659566063dSJacob Faibussowitsch   PetscCall(PetscPostIrecvInt(comm, tag2, nrqs, onodes2, olengths2, &rbuf2, &r_waits2));
266d5b485abSSatish Balay 
267d5b485abSSatish Balay   /*  Now  post the sends */
2689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrqr, &s_waits2));
269d5b485abSSatish Balay   for (i = 0; i < nrqr; ++i) {
270175c2bc5SJunchao Zhang     j = onodes1[i];
2719566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Isend(xdata[i], isz1[i], MPIU_INT, j, tag2, comm, s_waits2 + i));
272d5b485abSSatish Balay   }
273d5b485abSSatish Balay 
2749566063dSJacob Faibussowitsch   PetscCall(PetscFree(onodes1));
2759566063dSJacob Faibussowitsch   PetscCall(PetscFree(olengths1));
276175c2bc5SJunchao Zhang 
277d5b485abSSatish Balay   /* receive work done on other processors*/
278d5b485abSSatish Balay   {
279b24ad042SBarry Smith     PetscMPIInt idex;
280b24ad042SBarry Smith     PetscInt    is_no, ct1, max, *rbuf2_i, isz_i, *data_i, jmax;
281f1af5d2fSBarry Smith     PetscBT     table_i;
282d5b485abSSatish Balay 
283d5b485abSSatish Balay     for (i = 0; i < nrqs; ++i) {
2849566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Waitany(nrqs, r_waits2, &idex, MPI_STATUS_IGNORE));
285d5b485abSSatish Balay       /* Process the message*/
286999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
287d5b485abSSatish Balay       ct1     = 2 * rbuf2_i[0] + 1;
288999d9058SBarry Smith       jmax    = rbuf2[idex][0];
289d5b485abSSatish Balay       for (j = 1; j <= jmax; j++) {
290d5b485abSSatish Balay         max     = rbuf2_i[2 * j];
291d5b485abSSatish Balay         is_no   = rbuf2_i[2 * j - 1];
292d5b485abSSatish Balay         isz_i   = isz[is_no];
293d5b485abSSatish Balay         data_i  = data[is_no];
294d5b485abSSatish Balay         table_i = table[is_no];
295d5b485abSSatish Balay         for (k = 0; k < max; k++, ct1++) {
296d5b485abSSatish Balay           row = rbuf2_i[ct1];
29726fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row;
298d5b485abSSatish Balay         }
299d5b485abSSatish Balay         isz[is_no] = isz_i;
300d5b485abSSatish Balay       }
301d5b485abSSatish Balay     }
3029566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE));
303d5b485abSSatish Balay   }
304d5b485abSSatish Balay 
305d5b485abSSatish Balay   for (i = 0; i < imax; ++i) {
3069566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(iscomms[i], isz[i], data[i], PETSC_COPY_VALUES, is + i));
3079566063dSJacob Faibussowitsch     PetscCall(PetscCommDestroy(&iscomms[i]));
308d5b485abSSatish Balay   }
309d5b485abSSatish Balay 
3109566063dSJacob Faibussowitsch   PetscCall(PetscFree(iscomms));
3119566063dSJacob Faibussowitsch   PetscCall(PetscFree(onodes2));
3129566063dSJacob Faibussowitsch   PetscCall(PetscFree(olengths2));
313c7dd2924SSatish Balay 
3149566063dSJacob Faibussowitsch   PetscCall(PetscFree(pa));
315175c2bc5SJunchao Zhang   if (rbuf2) {
3169566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf2[0]));
3179566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf2));
318175c2bc5SJunchao Zhang   }
3199566063dSJacob Faibussowitsch   PetscCall(PetscFree(s_waits1));
3209566063dSJacob Faibussowitsch   PetscCall(PetscFree(r_waits1));
3219566063dSJacob Faibussowitsch   PetscCall(PetscFree(s_waits2));
3229566063dSJacob Faibussowitsch   PetscCall(PetscFree(r_waits2));
3239566063dSJacob Faibussowitsch   PetscCall(PetscFree5(table, data, isz, d_p, t_p));
324175c2bc5SJunchao Zhang   if (xdata) {
3259566063dSJacob Faibussowitsch     PetscCall(PetscFree(xdata[0]));
3269566063dSJacob Faibussowitsch     PetscCall(PetscFree(xdata));
327175c2bc5SJunchao Zhang   }
3289566063dSJacob Faibussowitsch   PetscCall(PetscFree(isz1));
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
330d5b485abSSatish Balay }
331d5b485abSSatish Balay 
332d5b485abSSatish Balay /*
333df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
334d5b485abSSatish Balay        the work on the local processor.
335d5b485abSSatish Balay 
336d5b485abSSatish Balay      Inputs:
337df36731bSBarry Smith       C      - MAT_MPIBAIJ;
338d5b485abSSatish Balay       imax - total no of index sets processed at a time;
339df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
340d5b485abSSatish Balay 
341d5b485abSSatish Balay      Output:
3420e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
343d5b485abSSatish Balay                to each index set;
344d5b485abSSatish Balay       data   - pointer to the solutions
345d5b485abSSatish Balay */
346d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C, PetscInt imax, PetscBT *table, PetscInt *isz, PetscInt **data)
347d71ae5a4SJacob Faibussowitsch {
348df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data;
349d5b485abSSatish Balay   Mat          A = c->A, B = c->B;
350df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
351b24ad042SBarry Smith   PetscInt     start, end, val, max, rstart, cstart, *ai, *aj;
352b24ad042SBarry Smith   PetscInt    *bi, *bj, *garray, i, j, k, row, *data_i, isz_i;
353f1af5d2fSBarry Smith   PetscBT      table_i;
354d5b485abSSatish Balay 
3553a40ed3dSBarry Smith   PetscFunctionBegin;
356899cda47SBarry Smith   rstart = c->rstartbs;
357899cda47SBarry Smith   cstart = c->cstartbs;
358d5b485abSSatish Balay   ai     = a->i;
359df36731bSBarry Smith   aj     = a->j;
360d5b485abSSatish Balay   bi     = b->i;
361df36731bSBarry Smith   bj     = b->j;
362d5b485abSSatish Balay   garray = c->garray;
363d5b485abSSatish Balay 
364d5b485abSSatish Balay   for (i = 0; i < imax; i++) {
365d5b485abSSatish Balay     data_i  = data[i];
366d5b485abSSatish Balay     table_i = table[i];
367d5b485abSSatish Balay     isz_i   = isz[i];
368d5b485abSSatish Balay     for (j = 0, max = isz[i]; j < max; j++) {
369d5b485abSSatish Balay       row   = data_i[j] - rstart;
370d5b485abSSatish Balay       start = ai[row];
371d5b485abSSatish Balay       end   = ai[row + 1];
372d5b485abSSatish Balay       for (k = start; k < end; k++) { /* Amat */
373df36731bSBarry Smith         val = aj[k] + cstart;
37426fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val;
375d5b485abSSatish Balay       }
376d5b485abSSatish Balay       start = bi[row];
377d5b485abSSatish Balay       end   = bi[row + 1];
378d5b485abSSatish Balay       for (k = start; k < end; k++) { /* Bmat */
379df36731bSBarry Smith         val = garray[bj[k]];
38026fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val;
381d5b485abSSatish Balay       }
382d5b485abSSatish Balay     }
383d5b485abSSatish Balay     isz[i] = isz_i;
384d5b485abSSatish Balay   }
3853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
386d5b485abSSatish Balay }
387d5b485abSSatish Balay /*
3882d4ee042Sprj-       MatIncreaseOverlap_MPIBAIJ_Receive - Process the received messages,
389d5b485abSSatish Balay          and return the output
390d5b485abSSatish Balay 
391d5b485abSSatish Balay          Input:
392d5b485abSSatish Balay            C    - the matrix
393d5b485abSSatish Balay            nrqr - no of messages being processed.
3942d4ee042Sprj-            rbuf - an array of pointers to the received requests
395d5b485abSSatish Balay 
396d5b485abSSatish Balay          Output:
397d5b485abSSatish Balay            xdata - array of messages to be sent back
398d5b485abSSatish Balay            isz1  - size of each message
399d5b485abSSatish Balay 
400a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
401d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
402a5b23f4aSJose E. Roman rather than all previous rows as it is now where a single large chunk of
403d5b485abSSatish Balay memory is used.
404d5b485abSSatish Balay 
405d5b485abSSatish Balay */
406d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C, PetscInt nrqr, PetscInt **rbuf, PetscInt **xdata, PetscInt *isz1)
407d71ae5a4SJacob Faibussowitsch {
408df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data;
409d5b485abSSatish Balay   Mat          A = c->A, B = c->B;
410df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
411b24ad042SBarry Smith   PetscInt     rstart, cstart, *ai, *aj, *bi, *bj, *garray, i, j, k;
412b24ad042SBarry Smith   PetscInt     row, total_sz, ct, ct1, ct2, ct3, mem_estimate, oct2, l, start, end;
413d2a212eaSBarry Smith   PetscInt     val, max1, max2, Mbs, no_malloc = 0, *tmp, new_estimate, ctr;
414b24ad042SBarry Smith   PetscInt    *rbuf_i, kmax, rbuf_0;
415f1af5d2fSBarry Smith   PetscBT      xtable;
416d5b485abSSatish Balay 
4173a40ed3dSBarry Smith   PetscFunctionBegin;
418df36731bSBarry Smith   Mbs    = c->Mbs;
419899cda47SBarry Smith   rstart = c->rstartbs;
420899cda47SBarry Smith   cstart = c->cstartbs;
421d5b485abSSatish Balay   ai     = a->i;
422df36731bSBarry Smith   aj     = a->j;
423d5b485abSSatish Balay   bi     = b->i;
424df36731bSBarry Smith   bj     = b->j;
425d5b485abSSatish Balay   garray = c->garray;
426d5b485abSSatish Balay 
427d5b485abSSatish Balay   for (i = 0, ct = 0, total_sz = 0; i < nrqr; ++i) {
428d5b485abSSatish Balay     rbuf_i = rbuf[i];
429d5b485abSSatish Balay     rbuf_0 = rbuf_i[0];
430d5b485abSSatish Balay     ct += rbuf_0;
43126fbe8dcSKarl Rupp     for (j = 1; j <= rbuf_0; j++) total_sz += rbuf_i[2 * j];
432d5b485abSSatish Balay   }
433d5b485abSSatish Balay 
434701b8814SBarry Smith   if (c->Mbs) max1 = ct * (a->nz + b->nz) / c->Mbs;
435701b8814SBarry Smith   else max1 = 1;
436d5b485abSSatish Balay   mem_estimate = 3 * ((total_sz > max1 ? total_sz : max1) + 1);
437175c2bc5SJunchao Zhang   if (nrqr) {
4389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mem_estimate, &xdata[0]));
439d5b485abSSatish Balay     ++no_malloc;
440175c2bc5SJunchao Zhang   }
4419566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(Mbs, &xtable));
4429566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(isz1, nrqr));
443d5b485abSSatish Balay 
444d5b485abSSatish Balay   ct3 = 0;
445d5b485abSSatish Balay   for (i = 0; i < nrqr; i++) { /* for easch mesg from proc i */
446d5b485abSSatish Balay     rbuf_i = rbuf[i];
447d5b485abSSatish Balay     rbuf_0 = rbuf_i[0];
448d5b485abSSatish Balay     ct1    = 2 * rbuf_0 + 1;
449d5b485abSSatish Balay     ct2    = ct1;
450d5b485abSSatish Balay     ct3 += ct1;
451d5b485abSSatish Balay     for (j = 1; j <= rbuf_0; j++) { /* for each IS from proc i*/
4529566063dSJacob Faibussowitsch       PetscCall(PetscBTMemzero(Mbs, xtable));
453d5b485abSSatish Balay       oct2 = ct2;
454d5b485abSSatish Balay       kmax = rbuf_i[2 * j];
455d5b485abSSatish Balay       for (k = 0; k < kmax; k++, ct1++) {
456d5b485abSSatish Balay         row = rbuf_i[ct1];
457f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable, row)) {
458d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
459b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5 * mem_estimate) + 1;
4609566063dSJacob Faibussowitsch             PetscCall(PetscMalloc1(new_estimate, &tmp));
4619566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate));
4629566063dSJacob Faibussowitsch             PetscCall(PetscFree(xdata[0]));
463d5b485abSSatish Balay             xdata[0]     = tmp;
4649371c9d4SSatish Balay             mem_estimate = new_estimate;
4659371c9d4SSatish Balay             ++no_malloc;
46626fbe8dcSKarl Rupp             for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1];
467d5b485abSSatish Balay           }
468d5b485abSSatish Balay           xdata[i][ct2++] = row;
469d5b485abSSatish Balay           ct3++;
470d5b485abSSatish Balay         }
471d5b485abSSatish Balay       }
472d5b485abSSatish Balay       for (k = oct2, max2 = ct2; k < max2; k++) {
473d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
474d5b485abSSatish Balay         start = ai[row];
475d5b485abSSatish Balay         end   = ai[row + 1];
476d5b485abSSatish Balay         for (l = start; l < end; l++) {
477df36731bSBarry Smith           val = aj[l] + cstart;
478f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable, val)) {
479d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
480b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5 * mem_estimate) + 1;
4819566063dSJacob Faibussowitsch               PetscCall(PetscMalloc1(new_estimate, &tmp));
4829566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate));
4839566063dSJacob Faibussowitsch               PetscCall(PetscFree(xdata[0]));
484d5b485abSSatish Balay               xdata[0]     = tmp;
4859371c9d4SSatish Balay               mem_estimate = new_estimate;
4869371c9d4SSatish Balay               ++no_malloc;
48726fbe8dcSKarl Rupp               for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1];
488d5b485abSSatish Balay             }
489d5b485abSSatish Balay             xdata[i][ct2++] = val;
490d5b485abSSatish Balay             ct3++;
491d5b485abSSatish Balay           }
492d5b485abSSatish Balay         }
493d5b485abSSatish Balay         start = bi[row];
494d5b485abSSatish Balay         end   = bi[row + 1];
495d5b485abSSatish Balay         for (l = start; l < end; l++) {
496df36731bSBarry Smith           val = garray[bj[l]];
497f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable, val)) {
498d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
499b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5 * mem_estimate) + 1;
5009566063dSJacob Faibussowitsch               PetscCall(PetscMalloc1(new_estimate, &tmp));
5019566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate));
5029566063dSJacob Faibussowitsch               PetscCall(PetscFree(xdata[0]));
503d5b485abSSatish Balay               xdata[0]     = tmp;
5049371c9d4SSatish Balay               mem_estimate = new_estimate;
5059371c9d4SSatish Balay               ++no_malloc;
50626fbe8dcSKarl Rupp               for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1];
507d5b485abSSatish Balay             }
508d5b485abSSatish Balay             xdata[i][ct2++] = val;
509d5b485abSSatish Balay             ct3++;
510d5b485abSSatish Balay           }
511d5b485abSSatish Balay         }
512d5b485abSSatish Balay       }
513d5b485abSSatish Balay       /* Update the header*/
514d5b485abSSatish Balay       xdata[i][2 * j]     = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
515d5b485abSSatish Balay       xdata[i][2 * j - 1] = rbuf_i[2 * j - 1];
516d5b485abSSatish Balay     }
517d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
518175c2bc5SJunchao Zhang     if (i + 1 < nrqr) xdata[i + 1] = xdata[i] + ct2;
519d5b485abSSatish Balay     isz1[i] = ct2; /* size of each message */
520d5b485abSSatish Balay   }
5219566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&xtable));
5229566063dSJacob Faibussowitsch   PetscCall(PetscInfo(C, "Allocated %" PetscInt_FMT " bytes, required %" PetscInt_FMT ", no of mallocs = %" PetscInt_FMT "\n", mem_estimate, ct3, no_malloc));
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
524d5b485abSSatish Balay }
525d5b485abSSatish Balay 
526d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submat[])
527d71ae5a4SJacob Faibussowitsch {
52817df9f7cSHong Zhang   IS          *isrow_block, *iscol_block;
529cf4f527aSSatish Balay   Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data;
530121971b2SHong Zhang   PetscInt     nmax, nstages_local, nstages, i, pos, max_no, N = C->cmap->N, bs = C->rmap->bs;
5315dc98d6aSHong Zhang   Mat_SeqBAIJ *subc;
5325c39f6d9SHong Zhang   Mat_SubSppt *smat;
533fdfbdca6SPierre Jolivet   PetscBool    sym = PETSC_FALSE, flg[2];
534a2feddc5SSatish Balay 
5353a40ed3dSBarry Smith   PetscFunctionBegin;
536fdfbdca6SPierre Jolivet   PetscCall(PetscObjectTypeCompare((PetscObject)C, MATMPISBAIJ, flg));
537fdfbdca6SPierre Jolivet   if (flg[0]) {
538fdfbdca6SPierre Jolivet     if (isrow == iscol) sym = PETSC_TRUE;
539fdfbdca6SPierre Jolivet     else {
540fdfbdca6SPierre Jolivet       flg[0] = flg[1] = PETSC_TRUE;
541fdfbdca6SPierre Jolivet       for (i = 0; i < ismax; i++) {
542fdfbdca6SPierre Jolivet         if (isrow[i] != iscol[i]) flg[0] = PETSC_FALSE;
543fdfbdca6SPierre Jolivet         PetscCall(ISGetLocalSize(iscol[0], &nmax));
544fdfbdca6SPierre Jolivet         if (nmax == C->cmap->N && flg[1]) PetscCall(ISIdentity(iscol[0], flg + 1));
545fdfbdca6SPierre Jolivet       }
546fdfbdca6SPierre Jolivet       sym = (PetscBool)(flg[0] || flg[1]);
547fdfbdca6SPierre Jolivet     }
548fdfbdca6SPierre Jolivet   }
54929dcf524SDmitry Karpeev   /* The compression and expansion should be avoided. Doesn't point
55029dcf524SDmitry Karpeev      out errors, might change the indices, hence buggey */
551864bd19bSPierre Jolivet   PetscCall(PetscMalloc2(ismax, &isrow_block, ismax, &iscol_block));
552864bd19bSPierre Jolivet   PetscCall(ISCompressIndicesGeneral(C->rmap->N, C->rmap->n, bs, ismax, isrow, isrow_block));
553864bd19bSPierre Jolivet   if (isrow == iscol) {
554864bd19bSPierre Jolivet     for (i = 0; i < ismax; i++) {
555864bd19bSPierre Jolivet       iscol_block[i] = isrow_block[i];
556864bd19bSPierre Jolivet       PetscCall(PetscObjectReference((PetscObject)iscol_block[i]));
557864bd19bSPierre Jolivet     }
558864bd19bSPierre Jolivet   } else PetscCall(ISCompressIndicesGeneral(N, C->cmap->n, bs, ismax, iscol, iscol_block));
559cf4f527aSSatish Balay 
560cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
5614698041cSHong Zhang   if (!C->cmap->N) nmax = 20 * 1000000 / sizeof(PetscInt);
5624698041cSHong Zhang   else nmax = 20 * 1000000 / (c->Nbs * sizeof(PetscInt));
563329f5518SBarry Smith   if (!nmax) nmax = 1;
5645dc98d6aSHong Zhang 
5655dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
5664698041cSHong Zhang     nstages_local = ismax / nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */
567cf4f527aSSatish Balay 
568653e4784SBarry Smith     /* Make sure every processor loops through the nstages */
5691c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&nstages_local, &nstages, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)C)));
5705dc98d6aSHong Zhang 
5714698041cSHong Zhang     /* Allocate memory to hold all the submatrices and dummy submatrices */
5729566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(ismax + nstages, submat));
5734698041cSHong Zhang   } else { /* MAT_REUSE_MATRIX */
5744698041cSHong Zhang     if (ismax) {
5755dc98d6aSHong Zhang       subc = (Mat_SeqBAIJ *)((*submat)[0]->data);
5765dc98d6aSHong Zhang       smat = subc->submatis1;
5774698041cSHong Zhang     } else { /* (*submat)[0] is a dummy matrix */
5785c39f6d9SHong Zhang       smat = (Mat_SubSppt *)(*submat)[0]->data;
5794698041cSHong Zhang     }
58028b400f6SJacob Faibussowitsch     PetscCheck(smat, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat");
5815dc98d6aSHong Zhang     nstages = smat->nstages;
5825dc98d6aSHong Zhang   }
5835dc98d6aSHong Zhang 
584cf4f527aSSatish Balay   for (i = 0, pos = 0; i < nstages; i++) {
585cf4f527aSSatish Balay     if (pos + nmax <= ismax) max_no = nmax;
5864e195584SJunchao Zhang     else if (pos >= ismax) max_no = 0;
587cf4f527aSSatish Balay     else max_no = ismax - pos;
588a42ab0d6SHong Zhang 
589fdfbdca6SPierre Jolivet     PetscCall(MatCreateSubMatrices_MPIBAIJ_local(C, max_no, isrow_block + pos, iscol_block + pos, scall, *submat + pos, sym));
5904e195584SJunchao Zhang     if (!max_no) {
5914e195584SJunchao Zhang       if (scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */
5925c39f6d9SHong Zhang         smat          = (Mat_SubSppt *)(*submat)[pos]->data;
5934698041cSHong Zhang         smat->nstages = nstages;
5944698041cSHong Zhang       }
5954e195584SJunchao Zhang       pos++; /* advance to next dummy matrix if any */
5964e195584SJunchao Zhang     } else pos += max_no;
597cf4f527aSSatish Balay   }
59836f4e84dSSatish Balay 
5995dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX && ismax) {
6005dc98d6aSHong Zhang     /* save nstages for reuse */
6015dc98d6aSHong Zhang     subc          = (Mat_SeqBAIJ *)((*submat)[0]->data);
6025dc98d6aSHong Zhang     smat          = subc->submatis1;
6035dc98d6aSHong Zhang     smat->nstages = nstages;
6045dc98d6aSHong Zhang   }
6055dc98d6aSHong Zhang 
60636f4e84dSSatish Balay   for (i = 0; i < ismax; i++) {
6079566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isrow_block[i]));
6089566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iscol_block[i]));
60936f4e84dSSatish Balay   }
6109566063dSJacob Faibussowitsch   PetscCall(PetscFree2(isrow_block, iscol_block));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
612a2feddc5SSatish Balay }
613a2feddc5SSatish Balay 
614b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
615fdfbdca6SPierre Jolivet PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submats, PetscBool sym)
616d71ae5a4SJacob Faibussowitsch {
61717df9f7cSHong Zhang   Mat_MPIBAIJ     *c = (Mat_MPIBAIJ *)C->data;
61817df9f7cSHong Zhang   Mat              A = c->A;
61917df9f7cSHong Zhang   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)c->B->data, *subc;
62017df9f7cSHong Zhang   const PetscInt **icol, **irow;
62117df9f7cSHong Zhang   PetscInt        *nrow, *ncol, start;
62217df9f7cSHong Zhang   PetscMPIInt      rank, size, tag0, tag2, tag3, tag4, *w1, *w2, *w3, *w4, nrqr;
623ddb3d6bcSHong Zhang   PetscInt       **sbuf1, **sbuf2, *sbuf2_i, i, j, k, l, ct1, ct2, **rbuf1, row, proc = -1;
62417df9f7cSHong Zhang   PetscInt         nrqs = 0, msz, **ptr = NULL, *req_size = NULL, *ctr = NULL, *pa, *tmp = NULL, tcol;
62517df9f7cSHong Zhang   PetscInt       **rbuf3 = NULL, *req_source1 = NULL, *req_source2, **sbuf_aj, **rbuf2 = NULL, max1, max2;
62617df9f7cSHong Zhang   PetscInt       **lens, is_no, ncols, *cols, mat_i, *mat_j, tmp2, jmax;
62717df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
628eec179cfSJacob Faibussowitsch   PetscHMapI *cmap, cmap_i = NULL, *rmap, rmap_i;
62917df9f7cSHong Zhang #else
63017df9f7cSHong Zhang   PetscInt **cmap, *cmap_i = NULL, **rmap, *rmap_i;
63117df9f7cSHong Zhang #endif
63296cff833SHong Zhang   const PetscInt *irow_i, *icol_i;
63317df9f7cSHong Zhang   PetscInt        ctr_j, *sbuf1_j, *sbuf_aj_i, *rbuf1_i, kmax, *lens_i;
63417df9f7cSHong Zhang   MPI_Request    *s_waits1, *r_waits1, *s_waits2, *r_waits2, *r_waits3;
63517df9f7cSHong Zhang   MPI_Request    *r_waits4, *s_waits3, *s_waits4;
63617df9f7cSHong Zhang   MPI_Comm        comm;
63764f5bb2dSSatish Balay   PetscScalar   **rbuf4, *rbuf4_i = NULL, **sbuf_aa, *vals, *mat_a = NULL, *imat_a = NULL, *sbuf_aa_i;
63817df9f7cSHong Zhang   PetscMPIInt    *onodes1, *olengths1, end;
63980d31651SHong Zhang   PetscInt      **row2proc, *row2proc_i, *imat_ilen, *imat_j, *imat_i;
6405c39f6d9SHong Zhang   Mat_SubSppt    *smat_i;
64117df9f7cSHong Zhang   PetscBool      *issorted, colflag, iscsorted = PETSC_TRUE;
64217df9f7cSHong Zhang   PetscInt       *sbuf1_i, *rbuf2_i, *rbuf3_i, ilen;
643a42ab0d6SHong Zhang   PetscInt        bs = C->rmap->bs, bs2 = c->bs2, rstart = c->rstartbs;
644a42ab0d6SHong Zhang   PetscBool       ijonly = c->ijonly; /* private flag indicates only matrix data structures are requested */
645a42ab0d6SHong Zhang   PetscInt        nzA, nzB, *a_i = a->i, *b_i = b->i, *a_j = a->j, *b_j = b->j, ctmp, imark, *cworkA, *cworkB;
646afb3cbd8SHong Zhang   PetscScalar    *vworkA = NULL, *vworkB = NULL, *a_a = a->a, *b_a = b->a;
647a42ab0d6SHong Zhang   PetscInt        cstart = c->cstartbs, *bmap = c->garray;
64880d31651SHong Zhang   PetscBool      *allrows, *allcolumns;
649a42ab0d6SHong Zhang 
65017df9f7cSHong Zhang   PetscFunctionBegin;
6519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
65217df9f7cSHong Zhang   size = c->size;
65317df9f7cSHong Zhang   rank = c->rank;
65417df9f7cSHong Zhang 
6559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(ismax, &row2proc, ismax, &cmap, ismax, &rmap, ismax + 1, &allcolumns, ismax, &allrows));
6569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(ismax, (PetscInt ***)&irow, ismax, (PetscInt ***)&icol, ismax, &nrow, ismax, &ncol, ismax, &issorted));
65717df9f7cSHong Zhang 
65817df9f7cSHong Zhang   for (i = 0; i < ismax; i++) {
6599566063dSJacob Faibussowitsch     PetscCall(ISSorted(iscol[i], &issorted[i]));
66080d31651SHong Zhang     if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */
6619566063dSJacob Faibussowitsch     PetscCall(ISSorted(isrow[i], &issorted[i]));
66217df9f7cSHong Zhang 
663ec3c739cSHong Zhang     /* Check for special case: allcolumns */
6649566063dSJacob Faibussowitsch     PetscCall(ISIdentity(iscol[i], &colflag));
6659566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(iscol[i], &ncol[i]));
6665dd0c08cSHong Zhang 
6671abc2f7bSHong Zhang     if (colflag && ncol[i] == c->Nbs) {
6681abc2f7bSHong Zhang       allcolumns[i] = PETSC_TRUE;
6691abc2f7bSHong Zhang       icol[i]       = NULL;
6701abc2f7bSHong Zhang     } else {
67117df9f7cSHong Zhang       allcolumns[i] = PETSC_FALSE;
6729566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(iscol[i], &icol[i]));
6731abc2f7bSHong Zhang     }
674ec3c739cSHong Zhang 
675ec3c739cSHong Zhang     /* Check for special case: allrows */
6769566063dSJacob Faibussowitsch     PetscCall(ISIdentity(isrow[i], &colflag));
6779566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(isrow[i], &nrow[i]));
678ec3c739cSHong Zhang     if (colflag && nrow[i] == c->Mbs) {
679ec3c739cSHong Zhang       allrows[i] = PETSC_TRUE;
680ec3c739cSHong Zhang       irow[i]    = NULL;
681ec3c739cSHong Zhang     } else {
682ec3c739cSHong Zhang       allrows[i] = PETSC_FALSE;
6839566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(isrow[i], &irow[i]));
684ec3c739cSHong Zhang     }
68517df9f7cSHong Zhang   }
68617df9f7cSHong Zhang 
68717df9f7cSHong Zhang   if (scall == MAT_REUSE_MATRIX) {
68817df9f7cSHong Zhang     /* Assumes new rows are same length as the old rows */
68917df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
690*f4f49eeaSPierre Jolivet       subc = (Mat_SeqBAIJ *)submats[i]->data;
691aed4548fSBarry Smith       PetscCheck(subc->mbs == nrow[i] && subc->nbs == ncol[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size");
69217df9f7cSHong Zhang 
69317df9f7cSHong Zhang       /* Initial matrix as if empty */
6949566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(subc->ilen, subc->mbs));
69517df9f7cSHong Zhang 
69617df9f7cSHong Zhang       /* Initial matrix as if empty */
69717df9f7cSHong Zhang       submats[i]->factortype = C->factortype;
69817df9f7cSHong Zhang 
69917df9f7cSHong Zhang       smat_i = subc->submatis1;
70017df9f7cSHong Zhang 
70117df9f7cSHong Zhang       nrqs        = smat_i->nrqs;
70217df9f7cSHong Zhang       nrqr        = smat_i->nrqr;
70317df9f7cSHong Zhang       rbuf1       = smat_i->rbuf1;
70417df9f7cSHong Zhang       rbuf2       = smat_i->rbuf2;
70517df9f7cSHong Zhang       rbuf3       = smat_i->rbuf3;
70617df9f7cSHong Zhang       req_source2 = smat_i->req_source2;
70717df9f7cSHong Zhang 
70817df9f7cSHong Zhang       sbuf1 = smat_i->sbuf1;
70917df9f7cSHong Zhang       sbuf2 = smat_i->sbuf2;
71017df9f7cSHong Zhang       ptr   = smat_i->ptr;
71117df9f7cSHong Zhang       tmp   = smat_i->tmp;
71217df9f7cSHong Zhang       ctr   = smat_i->ctr;
71317df9f7cSHong Zhang 
71417df9f7cSHong Zhang       pa          = smat_i->pa;
71517df9f7cSHong Zhang       req_size    = smat_i->req_size;
71617df9f7cSHong Zhang       req_source1 = smat_i->req_source1;
71717df9f7cSHong Zhang 
71817df9f7cSHong Zhang       allcolumns[i] = smat_i->allcolumns;
719ec3c739cSHong Zhang       allrows[i]    = smat_i->allrows;
72017df9f7cSHong Zhang       row2proc[i]   = smat_i->row2proc;
72117df9f7cSHong Zhang       rmap[i]       = smat_i->rmap;
72217df9f7cSHong Zhang       cmap[i]       = smat_i->cmap;
72317df9f7cSHong Zhang     }
7244698041cSHong Zhang 
7254698041cSHong Zhang     if (!ismax) { /* Get dummy submatrices and retrieve struct submatis1 */
72608401ef6SPierre Jolivet       PetscCheck(submats[0], PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "submats are null, cannot reuse");
7275c39f6d9SHong Zhang       smat_i = (Mat_SubSppt *)submats[0]->data;
7284698041cSHong Zhang 
7294698041cSHong Zhang       nrqs        = smat_i->nrqs;
7304698041cSHong Zhang       nrqr        = smat_i->nrqr;
7314698041cSHong Zhang       rbuf1       = smat_i->rbuf1;
7324698041cSHong Zhang       rbuf2       = smat_i->rbuf2;
7334698041cSHong Zhang       rbuf3       = smat_i->rbuf3;
7344698041cSHong Zhang       req_source2 = smat_i->req_source2;
7354698041cSHong Zhang 
7364698041cSHong Zhang       sbuf1 = smat_i->sbuf1;
7374698041cSHong Zhang       sbuf2 = smat_i->sbuf2;
7384698041cSHong Zhang       ptr   = smat_i->ptr;
7394698041cSHong Zhang       tmp   = smat_i->tmp;
7404698041cSHong Zhang       ctr   = smat_i->ctr;
7414698041cSHong Zhang 
7424698041cSHong Zhang       pa          = smat_i->pa;
7434698041cSHong Zhang       req_size    = smat_i->req_size;
7444698041cSHong Zhang       req_source1 = smat_i->req_source1;
7454698041cSHong Zhang 
746c67fe082SHong Zhang       allcolumns[0] = PETSC_FALSE;
7474698041cSHong Zhang     }
74817df9f7cSHong Zhang   } else { /* scall == MAT_INITIAL_MATRIX */
74917df9f7cSHong Zhang     /* Get some new tags to keep the communication clean */
7509566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2));
7519566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag3));
75217df9f7cSHong Zhang 
75317df9f7cSHong Zhang     /* evaluate communication - mesg to who, length of mesg, and buffer space
75417df9f7cSHong Zhang      required. Based on this, buffers are allocated, and data copied into them*/
7559566063dSJacob Faibussowitsch     PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4)); /* mesg size, initialize work vectors */
75617df9f7cSHong Zhang 
75717df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
75817df9f7cSHong Zhang       jmax   = nrow[i];
75917df9f7cSHong Zhang       irow_i = irow[i];
76017df9f7cSHong Zhang 
7619566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jmax, &row2proc_i));
76217df9f7cSHong Zhang       row2proc[i] = row2proc_i;
76317df9f7cSHong Zhang 
76417df9f7cSHong Zhang       if (issorted[i]) proc = 0;
76517df9f7cSHong Zhang       for (j = 0; j < jmax; j++) {
76617df9f7cSHong Zhang         if (!issorted[i]) proc = 0;
767ec3c739cSHong Zhang         if (allrows[i]) row = j;
768ec3c739cSHong Zhang         else row = irow_i[j];
769ec3c739cSHong Zhang 
770a42ab0d6SHong Zhang         while (row >= c->rangebs[proc + 1]) proc++;
77117df9f7cSHong Zhang         w4[proc]++;
77217df9f7cSHong Zhang         row2proc_i[j] = proc; /* map row index to proc */
77317df9f7cSHong Zhang       }
77417df9f7cSHong Zhang       for (j = 0; j < size; j++) {
7759371c9d4SSatish Balay         if (w4[j]) {
7769371c9d4SSatish Balay           w1[j] += w4[j];
7779371c9d4SSatish Balay           w3[j]++;
7789371c9d4SSatish Balay           w4[j] = 0;
7799371c9d4SSatish Balay         }
78017df9f7cSHong Zhang       }
78117df9f7cSHong Zhang     }
78217df9f7cSHong Zhang 
78317df9f7cSHong Zhang     nrqs     = 0; /* no of outgoing messages */
78417df9f7cSHong Zhang     msz      = 0; /* total mesg length (for all procs) */
78517df9f7cSHong Zhang     w1[rank] = 0; /* no mesg sent to self */
78617df9f7cSHong Zhang     w3[rank] = 0;
78717df9f7cSHong Zhang     for (i = 0; i < size; i++) {
7889371c9d4SSatish Balay       if (w1[i]) {
7899371c9d4SSatish Balay         w2[i] = 1;
7909371c9d4SSatish Balay         nrqs++;
7919371c9d4SSatish Balay       } /* there exists a message to proc i */
79217df9f7cSHong Zhang     }
7939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &pa)); /*(proc -array)*/
79417df9f7cSHong Zhang     for (i = 0, j = 0; i < size; i++) {
7959371c9d4SSatish Balay       if (w1[i]) {
7969371c9d4SSatish Balay         pa[j] = i;
7979371c9d4SSatish Balay         j++;
7989371c9d4SSatish Balay       }
79917df9f7cSHong Zhang     }
80017df9f7cSHong Zhang 
80117df9f7cSHong Zhang     /* Each message would have a header = 1 + 2*(no of IS) + data */
80217df9f7cSHong Zhang     for (i = 0; i < nrqs; i++) {
80317df9f7cSHong Zhang       j = pa[i];
80417df9f7cSHong Zhang       w1[j] += w2[j] + 2 * w3[j];
80517df9f7cSHong Zhang       msz += w1[j];
80617df9f7cSHong Zhang     }
8079566063dSJacob Faibussowitsch     PetscCall(PetscInfo(0, "Number of outgoing messages %" PetscInt_FMT " Total message length %" PetscInt_FMT "\n", nrqs, msz));
80817df9f7cSHong Zhang 
80917df9f7cSHong Zhang     /* Determine the number of messages to expect, their lengths, from from-ids */
8109566063dSJacob Faibussowitsch     PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr));
8119566063dSJacob Faibussowitsch     PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1));
81217df9f7cSHong Zhang 
81317df9f7cSHong Zhang     /* Now post the Irecvs corresponding to these messages */
8149566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag0));
8159566063dSJacob Faibussowitsch     PetscCall(PetscPostIrecvInt(comm, tag0, nrqr, onodes1, olengths1, &rbuf1, &r_waits1));
81617df9f7cSHong Zhang 
81717df9f7cSHong Zhang     /* Allocate Memory for outgoing messages */
8189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(size, &sbuf1, size, &ptr, 2 * msz, &tmp, size, &ctr));
8199566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(sbuf1, size));
8209566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ptr, size));
82117df9f7cSHong Zhang 
82217df9f7cSHong Zhang     {
82317df9f7cSHong Zhang       PetscInt *iptr = tmp;
82417df9f7cSHong Zhang       k              = 0;
82517df9f7cSHong Zhang       for (i = 0; i < nrqs; i++) {
82617df9f7cSHong Zhang         j = pa[i];
82717df9f7cSHong Zhang         iptr += k;
82817df9f7cSHong Zhang         sbuf1[j] = iptr;
82917df9f7cSHong Zhang         k        = w1[j];
83017df9f7cSHong Zhang       }
83117df9f7cSHong Zhang     }
83217df9f7cSHong Zhang 
83317df9f7cSHong Zhang     /* Form the outgoing messages. Initialize the header space */
83417df9f7cSHong Zhang     for (i = 0; i < nrqs; i++) {
83517df9f7cSHong Zhang       j           = pa[i];
83617df9f7cSHong Zhang       sbuf1[j][0] = 0;
8379566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(sbuf1[j] + 1, 2 * w3[j]));
83817df9f7cSHong Zhang       ptr[j] = sbuf1[j] + 2 * w3[j] + 1;
83917df9f7cSHong Zhang     }
84017df9f7cSHong Zhang 
84117df9f7cSHong Zhang     /* Parse the isrow and copy data into outbuf */
84217df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
84317df9f7cSHong Zhang       row2proc_i = row2proc[i];
8449566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctr, size));
84517df9f7cSHong Zhang       irow_i = irow[i];
84617df9f7cSHong Zhang       jmax   = nrow[i];
84717df9f7cSHong Zhang       for (j = 0; j < jmax; j++) { /* parse the indices of each IS */
84817df9f7cSHong Zhang         proc = row2proc_i[j];
849ec3c739cSHong Zhang         if (allrows[i]) row = j;
850ec3c739cSHong Zhang         else row = irow_i[j];
851ec3c739cSHong Zhang 
85217df9f7cSHong Zhang         if (proc != rank) { /* copy to the outgoing buf*/
85317df9f7cSHong Zhang           ctr[proc]++;
854ec3c739cSHong Zhang           *ptr[proc] = row;
85517df9f7cSHong Zhang           ptr[proc]++;
85617df9f7cSHong Zhang         }
85717df9f7cSHong Zhang       }
85817df9f7cSHong Zhang       /* Update the headers for the current IS */
85917df9f7cSHong Zhang       for (j = 0; j < size; j++) { /* Can Optimise this loop too */
86017df9f7cSHong Zhang         if ((ctr_j = ctr[j])) {
86117df9f7cSHong Zhang           sbuf1_j            = sbuf1[j];
86217df9f7cSHong Zhang           k                  = ++sbuf1_j[0];
86317df9f7cSHong Zhang           sbuf1_j[2 * k]     = ctr_j;
86417df9f7cSHong Zhang           sbuf1_j[2 * k - 1] = i;
86517df9f7cSHong Zhang         }
86617df9f7cSHong Zhang       }
86717df9f7cSHong Zhang     }
86817df9f7cSHong Zhang 
86917df9f7cSHong Zhang     /*  Now  post the sends */
8709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &s_waits1));
87117df9f7cSHong Zhang     for (i = 0; i < nrqs; ++i) {
87217df9f7cSHong Zhang       j = pa[i];
8739566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sbuf1[j], w1[j], MPIU_INT, j, tag0, comm, s_waits1 + i));
87417df9f7cSHong Zhang     }
87517df9f7cSHong Zhang 
87617df9f7cSHong Zhang     /* Post Receives to capture the buffer size */
8779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &r_waits2));
8789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(nrqs, &req_source2, nrqs, &rbuf2, nrqs, &rbuf3));
879175c2bc5SJunchao Zhang     if (nrqs) rbuf2[0] = tmp + msz;
880ad540459SPierre Jolivet     for (i = 1; i < nrqs; ++i) rbuf2[i] = rbuf2[i - 1] + w1[pa[i - 1]];
88117df9f7cSHong Zhang     for (i = 0; i < nrqs; ++i) {
88217df9f7cSHong Zhang       j = pa[i];
8839566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(rbuf2[i], w1[j], MPIU_INT, j, tag2, comm, r_waits2 + i));
88417df9f7cSHong Zhang     }
88517df9f7cSHong Zhang 
88617df9f7cSHong Zhang     /* Send to other procs the buf size they should allocate */
88717df9f7cSHong Zhang     /* Receive messages*/
8889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &s_waits2));
8899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(nrqr, &sbuf2, nrqr, &req_size, nrqr, &req_source1));
89017df9f7cSHong Zhang 
8919566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE));
89217df9f7cSHong Zhang     for (i = 0; i < nrqr; ++i) {
89317df9f7cSHong Zhang       req_size[i] = 0;
89417df9f7cSHong Zhang       rbuf1_i     = rbuf1[i];
89517df9f7cSHong Zhang       start       = 2 * rbuf1_i[0] + 1;
896175c2bc5SJunchao Zhang       end         = olengths1[i];
8979566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(end, &sbuf2[i]));
89817df9f7cSHong Zhang       sbuf2_i = sbuf2[i];
89917df9f7cSHong Zhang       for (j = start; j < end; j++) {
900ddb3d6bcSHong Zhang         row        = rbuf1_i[j] - rstart;
901ddb3d6bcSHong Zhang         ncols      = a_i[row + 1] - a_i[row] + b_i[row + 1] - b_i[row];
90217df9f7cSHong Zhang         sbuf2_i[j] = ncols;
90317df9f7cSHong Zhang         req_size[i] += ncols;
90417df9f7cSHong Zhang       }
905175c2bc5SJunchao Zhang       req_source1[i] = onodes1[i];
90617df9f7cSHong Zhang       /* form the header */
90717df9f7cSHong Zhang       sbuf2_i[0] = req_size[i];
90817df9f7cSHong Zhang       for (j = 1; j < start; j++) sbuf2_i[j] = rbuf1_i[j];
90917df9f7cSHong Zhang 
9109566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sbuf2_i, end, MPIU_INT, req_source1[i], tag2, comm, s_waits2 + i));
91117df9f7cSHong Zhang     }
912ddb3d6bcSHong Zhang 
9139566063dSJacob Faibussowitsch     PetscCall(PetscFree(onodes1));
9149566063dSJacob Faibussowitsch     PetscCall(PetscFree(olengths1));
915175c2bc5SJunchao Zhang 
9169566063dSJacob Faibussowitsch     PetscCall(PetscFree(r_waits1));
9179566063dSJacob Faibussowitsch     PetscCall(PetscFree4(w1, w2, w3, w4));
91817df9f7cSHong Zhang 
91917df9f7cSHong Zhang     /* Receive messages*/
9209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &r_waits3));
92117df9f7cSHong Zhang 
9229566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqs, r_waits2, MPI_STATUSES_IGNORE));
92317df9f7cSHong Zhang     for (i = 0; i < nrqs; ++i) {
9249566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(rbuf2[i][0], &rbuf3[i]));
925175c2bc5SJunchao Zhang       req_source2[i] = pa[i];
9269566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(rbuf3[i], rbuf2[i][0], MPIU_INT, req_source2[i], tag3, comm, r_waits3 + i));
92717df9f7cSHong Zhang     }
9289566063dSJacob Faibussowitsch     PetscCall(PetscFree(r_waits2));
92917df9f7cSHong Zhang 
93017df9f7cSHong Zhang     /* Wait on sends1 and sends2 */
9319566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE));
9329566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE));
9339566063dSJacob Faibussowitsch     PetscCall(PetscFree(s_waits1));
9349566063dSJacob Faibussowitsch     PetscCall(PetscFree(s_waits2));
93517df9f7cSHong Zhang 
93617df9f7cSHong Zhang     /* Now allocate sending buffers for a->j, and send them off */
9379566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &sbuf_aj));
93817df9f7cSHong Zhang     for (i = 0, j = 0; i < nrqr; i++) j += req_size[i];
9399566063dSJacob Faibussowitsch     if (nrqr) PetscCall(PetscMalloc1(j, &sbuf_aj[0]));
94017df9f7cSHong Zhang     for (i = 1; i < nrqr; i++) sbuf_aj[i] = sbuf_aj[i - 1] + req_size[i - 1];
94117df9f7cSHong Zhang 
9429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &s_waits3));
94317df9f7cSHong Zhang     {
94417df9f7cSHong Zhang       for (i = 0; i < nrqr; i++) {
94517df9f7cSHong Zhang         rbuf1_i   = rbuf1[i];
94617df9f7cSHong Zhang         sbuf_aj_i = sbuf_aj[i];
94717df9f7cSHong Zhang         ct1       = 2 * rbuf1_i[0] + 1;
94817df9f7cSHong Zhang         ct2       = 0;
94917df9f7cSHong Zhang         for (j = 1, max1 = rbuf1_i[0]; j <= max1; j++) {
95017df9f7cSHong Zhang           kmax = rbuf1[i][2 * j];
95117df9f7cSHong Zhang           for (k = 0; k < kmax; k++, ct1++) {
95217df9f7cSHong Zhang             row    = rbuf1_i[ct1] - rstart;
9539371c9d4SSatish Balay             nzA    = a_i[row + 1] - a_i[row];
9549371c9d4SSatish Balay             nzB    = b_i[row + 1] - b_i[row];
95517df9f7cSHong Zhang             ncols  = nzA + nzB;
9568e3a54c0SPierre Jolivet             cworkA = PetscSafePointerPlusOffset(a_j, a_i[row]);
9578e3a54c0SPierre Jolivet             cworkB = PetscSafePointerPlusOffset(b_j, b_i[row]);
95817df9f7cSHong Zhang 
95917df9f7cSHong Zhang             /* load the column indices for this row into cols */
96017df9f7cSHong Zhang             cols = sbuf_aj_i + ct2;
96117df9f7cSHong Zhang             for (l = 0; l < nzB; l++) {
962a42ab0d6SHong Zhang               if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
963a42ab0d6SHong Zhang               else break;
96417df9f7cSHong Zhang             }
965a42ab0d6SHong Zhang             imark = l;
966ad540459SPierre Jolivet             for (l = 0; l < nzA; l++) cols[imark + l] = cstart + cworkA[l];
967a42ab0d6SHong Zhang             for (l = imark; l < nzB; l++) cols[nzA + l] = bmap[cworkB[l]];
96817df9f7cSHong Zhang             ct2 += ncols;
96917df9f7cSHong Zhang           }
97017df9f7cSHong Zhang         }
9719566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Isend(sbuf_aj_i, req_size[i], MPIU_INT, req_source1[i], tag3, comm, s_waits3 + i));
97217df9f7cSHong Zhang       }
97317df9f7cSHong Zhang     }
97417df9f7cSHong Zhang 
97517df9f7cSHong Zhang     /* create col map: global col of C -> local col of submatrices */
97617df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
97717df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
97817df9f7cSHong Zhang       if (!allcolumns[i]) {
979eec179cfSJacob Faibussowitsch         PetscCall(PetscHMapICreateWithSize(ncol[i], cmap + i));
98017df9f7cSHong Zhang 
98117df9f7cSHong Zhang         jmax   = ncol[i];
98217df9f7cSHong Zhang         icol_i = icol[i];
98317df9f7cSHong Zhang         cmap_i = cmap[i];
984c76ffc5fSJacob Faibussowitsch         for (j = 0; j < jmax; j++) PetscCall(PetscHMapISet(cmap[i], icol_i[j] + 1, j + 1));
98517df9f7cSHong Zhang       } else cmap[i] = NULL;
98617df9f7cSHong Zhang     }
98717df9f7cSHong Zhang #else
98817df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
98917df9f7cSHong Zhang       if (!allcolumns[i]) {
9909566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(c->Nbs, &cmap[i]));
99117df9f7cSHong Zhang         jmax   = ncol[i];
99217df9f7cSHong Zhang         icol_i = icol[i];
99317df9f7cSHong Zhang         cmap_i = cmap[i];
994a42ab0d6SHong Zhang         for (j = 0; j < jmax; j++) cmap_i[icol_i[j]] = j + 1;
99517df9f7cSHong Zhang       } else cmap[i] = NULL;
99617df9f7cSHong Zhang     }
99717df9f7cSHong Zhang #endif
99817df9f7cSHong Zhang 
99917df9f7cSHong Zhang     /* Create lens which is required for MatCreate... */
100017df9f7cSHong Zhang     for (i = 0, j = 0; i < ismax; i++) j += nrow[i];
10019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ismax, &lens));
100217df9f7cSHong Zhang 
100348a46eb9SPierre Jolivet     if (ismax) PetscCall(PetscCalloc1(j, &lens[0]));
10048e3a54c0SPierre Jolivet     for (i = 1; i < ismax; i++) lens[i] = PetscSafePointerPlusOffset(lens[i - 1], nrow[i - 1]);
100517df9f7cSHong Zhang 
100617df9f7cSHong Zhang     /* Update lens from local data */
100717df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
100817df9f7cSHong Zhang       row2proc_i = row2proc[i];
100917df9f7cSHong Zhang       jmax       = nrow[i];
101017df9f7cSHong Zhang       if (!allcolumns[i]) cmap_i = cmap[i];
101117df9f7cSHong Zhang       irow_i = irow[i];
101217df9f7cSHong Zhang       lens_i = lens[i];
101317df9f7cSHong Zhang       for (j = 0; j < jmax; j++) {
1014ec3c739cSHong Zhang         if (allrows[i]) row = j;
1015ec3c739cSHong Zhang         else row = irow_i[j]; /* global blocked row of C */
1016ec3c739cSHong Zhang 
101717df9f7cSHong Zhang         proc = row2proc_i[j];
101817df9f7cSHong Zhang         if (proc == rank) {
1019a42ab0d6SHong Zhang           /* Get indices from matA and then from matB */
102017df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1021a42ab0d6SHong Zhang           PetscInt tt;
1022a42ab0d6SHong Zhang #endif
1023a42ab0d6SHong Zhang           row    = row - rstart;
1024a42ab0d6SHong Zhang           nzA    = a_i[row + 1] - a_i[row];
1025a42ab0d6SHong Zhang           nzB    = b_i[row + 1] - b_i[row];
1026a42ab0d6SHong Zhang           cworkA = a_j + a_i[row];
10278e3a54c0SPierre Jolivet           cworkB = PetscSafePointerPlusOffset(b_j, b_i[row]);
1028a42ab0d6SHong Zhang 
1029a42ab0d6SHong Zhang           if (!allcolumns[i]) {
1030a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1031a42ab0d6SHong Zhang             for (k = 0; k < nzA; k++) {
1032eec179cfSJacob Faibussowitsch               PetscCall(PetscHMapIGetWithDefault(cmap_i, cstart + cworkA[k] + 1, 0, &tt));
1033a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1034a42ab0d6SHong Zhang             }
1035a42ab0d6SHong Zhang             for (k = 0; k < nzB; k++) {
1036eec179cfSJacob Faibussowitsch               PetscCall(PetscHMapIGetWithDefault(cmap_i, bmap[cworkB[k]] + 1, 0, &tt));
1037a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1038a42ab0d6SHong Zhang             }
1039a42ab0d6SHong Zhang 
1040a42ab0d6SHong Zhang #else
1041a42ab0d6SHong Zhang             for (k = 0; k < nzA; k++) {
1042a42ab0d6SHong Zhang               if (cmap_i[cstart + cworkA[k]]) lens_i[j]++;
1043a42ab0d6SHong Zhang             }
1044a42ab0d6SHong Zhang             for (k = 0; k < nzB; k++) {
1045a42ab0d6SHong Zhang               if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++;
1046a42ab0d6SHong Zhang             }
1047a42ab0d6SHong Zhang #endif
1048a42ab0d6SHong Zhang           } else { /* allcolumns */
1049a42ab0d6SHong Zhang             lens_i[j] = nzA + nzB;
1050a42ab0d6SHong Zhang           }
105117df9f7cSHong Zhang         }
105217df9f7cSHong Zhang       }
105317df9f7cSHong Zhang     }
105417df9f7cSHong Zhang 
105517df9f7cSHong Zhang     /* Create row map: global row of C -> local row of submatrices */
105617df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
1057059f6b74SHong Zhang       if (!allrows[i]) {
1058059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE)
1059eec179cfSJacob Faibussowitsch         PetscCall(PetscHMapICreateWithSize(nrow[i], rmap + i));
106017df9f7cSHong Zhang         irow_i = irow[i];
106117df9f7cSHong Zhang         jmax   = nrow[i];
106217df9f7cSHong Zhang         for (j = 0; j < jmax; j++) {
1063ec3c739cSHong Zhang           if (allrows[i]) {
1064c76ffc5fSJacob Faibussowitsch             PetscCall(PetscHMapISet(rmap[i], j + 1, j + 1));
1065ec3c739cSHong Zhang           } else {
1066c76ffc5fSJacob Faibussowitsch             PetscCall(PetscHMapISet(rmap[i], irow_i[j] + 1, j + 1));
106717df9f7cSHong Zhang           }
106817df9f7cSHong Zhang         }
106917df9f7cSHong Zhang #else
10709566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(c->Mbs, &rmap[i]));
107117df9f7cSHong Zhang         rmap_i = rmap[i];
107217df9f7cSHong Zhang         irow_i = irow[i];
107317df9f7cSHong Zhang         jmax   = nrow[i];
107417df9f7cSHong Zhang         for (j = 0; j < jmax; j++) {
1075ec3c739cSHong Zhang           if (allrows[i]) rmap_i[j] = j;
1076ec3c739cSHong Zhang           else rmap_i[irow_i[j]] = j;
107717df9f7cSHong Zhang         }
107817df9f7cSHong Zhang #endif
1079059f6b74SHong Zhang       } else rmap[i] = NULL;
1080059f6b74SHong Zhang     }
108117df9f7cSHong Zhang 
108217df9f7cSHong Zhang     /* Update lens from offproc data */
108317df9f7cSHong Zhang     {
108417df9f7cSHong Zhang       PetscInt *rbuf2_i, *rbuf3_i, *sbuf1_i;
108517df9f7cSHong Zhang 
10869566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Waitall(nrqs, r_waits3, MPI_STATUSES_IGNORE));
108717df9f7cSHong Zhang       for (tmp2 = 0; tmp2 < nrqs; tmp2++) {
108817df9f7cSHong Zhang         sbuf1_i = sbuf1[pa[tmp2]];
108917df9f7cSHong Zhang         jmax    = sbuf1_i[0];
109017df9f7cSHong Zhang         ct1     = 2 * jmax + 1;
109117df9f7cSHong Zhang         ct2     = 0;
109217df9f7cSHong Zhang         rbuf2_i = rbuf2[tmp2];
109317df9f7cSHong Zhang         rbuf3_i = rbuf3[tmp2];
109417df9f7cSHong Zhang         for (j = 1; j <= jmax; j++) {
109517df9f7cSHong Zhang           is_no  = sbuf1_i[2 * j - 1];
109617df9f7cSHong Zhang           max1   = sbuf1_i[2 * j];
109717df9f7cSHong Zhang           lens_i = lens[is_no];
109817df9f7cSHong Zhang           if (!allcolumns[is_no]) cmap_i = cmap[is_no];
109917df9f7cSHong Zhang           rmap_i = rmap[is_no];
110017df9f7cSHong Zhang           for (k = 0; k < max1; k++, ct1++) {
1101059f6b74SHong Zhang             if (allrows[is_no]) {
1102059f6b74SHong Zhang               row = sbuf1_i[ct1];
1103059f6b74SHong Zhang             } else {
110417df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1105eec179cfSJacob Faibussowitsch               PetscCall(PetscHMapIGetWithDefault(rmap_i, sbuf1_i[ct1] + 1, 0, &row));
110617df9f7cSHong Zhang               row--;
110708401ef6SPierre Jolivet               PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table");
110817df9f7cSHong Zhang #else
110917df9f7cSHong Zhang               row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
111017df9f7cSHong Zhang #endif
1111059f6b74SHong Zhang             }
111217df9f7cSHong Zhang             max2 = rbuf2_i[ct1];
111317df9f7cSHong Zhang             for (l = 0; l < max2; l++, ct2++) {
111417df9f7cSHong Zhang               if (!allcolumns[is_no]) {
111517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1116eec179cfSJacob Faibussowitsch                 PetscCall(PetscHMapIGetWithDefault(cmap_i, rbuf3_i[ct2] + 1, 0, &tcol));
111717df9f7cSHong Zhang #else
111817df9f7cSHong Zhang                 tcol = cmap_i[rbuf3_i[ct2]];
111917df9f7cSHong Zhang #endif
112017df9f7cSHong Zhang                 if (tcol) lens_i[row]++;
112117df9f7cSHong Zhang               } else {         /* allcolumns */
112217df9f7cSHong Zhang                 lens_i[row]++; /* lens_i[row] += max2 ? */
112317df9f7cSHong Zhang               }
112417df9f7cSHong Zhang             }
112517df9f7cSHong Zhang           }
112617df9f7cSHong Zhang         }
112717df9f7cSHong Zhang       }
112817df9f7cSHong Zhang     }
11299566063dSJacob Faibussowitsch     PetscCall(PetscFree(r_waits3));
11309566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, s_waits3, MPI_STATUSES_IGNORE));
11319566063dSJacob Faibussowitsch     PetscCall(PetscFree(s_waits3));
113217df9f7cSHong Zhang 
113317df9f7cSHong Zhang     /* Create the submatrices */
113417df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
1135a42ab0d6SHong Zhang       PetscInt bs_tmp;
1136a42ab0d6SHong Zhang       if (ijonly) bs_tmp = 1;
1137a42ab0d6SHong Zhang       else bs_tmp = bs;
113817df9f7cSHong Zhang 
11399566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, submats + i));
11409566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(submats[i], nrow[i] * bs_tmp, ncol[i] * bs_tmp, PETSC_DETERMINE, PETSC_DETERMINE));
114117df9f7cSHong Zhang 
1142fdfbdca6SPierre Jolivet       PetscCall(MatSetType(submats[i], sym ? ((PetscObject)A)->type_name : MATSEQBAIJ));
11439566063dSJacob Faibussowitsch       PetscCall(MatSeqBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i]));
11449566063dSJacob Faibussowitsch       PetscCall(MatSeqSBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i])); /* this subroutine is used by SBAIJ routines */
114517df9f7cSHong Zhang 
11465c39f6d9SHong Zhang       /* create struct Mat_SubSppt and attached it to submat */
11479566063dSJacob Faibussowitsch       PetscCall(PetscNew(&smat_i));
114817df9f7cSHong Zhang       subc            = (Mat_SeqBAIJ *)submats[i]->data;
114917df9f7cSHong Zhang       subc->submatis1 = smat_i;
115017df9f7cSHong Zhang 
115117df9f7cSHong Zhang       smat_i->destroy          = submats[i]->ops->destroy;
1152f68bb481SHong Zhang       submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ;
115317df9f7cSHong Zhang       submats[i]->factortype   = C->factortype;
115417df9f7cSHong Zhang 
115517df9f7cSHong Zhang       smat_i->id          = i;
115617df9f7cSHong Zhang       smat_i->nrqs        = nrqs;
115717df9f7cSHong Zhang       smat_i->nrqr        = nrqr;
115817df9f7cSHong Zhang       smat_i->rbuf1       = rbuf1;
115917df9f7cSHong Zhang       smat_i->rbuf2       = rbuf2;
116017df9f7cSHong Zhang       smat_i->rbuf3       = rbuf3;
116117df9f7cSHong Zhang       smat_i->sbuf2       = sbuf2;
116217df9f7cSHong Zhang       smat_i->req_source2 = req_source2;
116317df9f7cSHong Zhang 
116417df9f7cSHong Zhang       smat_i->sbuf1 = sbuf1;
116517df9f7cSHong Zhang       smat_i->ptr   = ptr;
116617df9f7cSHong Zhang       smat_i->tmp   = tmp;
116717df9f7cSHong Zhang       smat_i->ctr   = ctr;
116817df9f7cSHong Zhang 
116917df9f7cSHong Zhang       smat_i->pa          = pa;
117017df9f7cSHong Zhang       smat_i->req_size    = req_size;
117117df9f7cSHong Zhang       smat_i->req_source1 = req_source1;
117217df9f7cSHong Zhang 
117317df9f7cSHong Zhang       smat_i->allcolumns = allcolumns[i];
1174ec3c739cSHong Zhang       smat_i->allrows    = allrows[i];
117517df9f7cSHong Zhang       smat_i->singleis   = PETSC_FALSE;
117617df9f7cSHong Zhang       smat_i->row2proc   = row2proc[i];
117717df9f7cSHong Zhang       smat_i->rmap       = rmap[i];
117817df9f7cSHong Zhang       smat_i->cmap       = cmap[i];
117917df9f7cSHong Zhang     }
118017df9f7cSHong Zhang 
11814698041cSHong Zhang     if (!ismax) { /* Create dummy submats[0] for reuse struct subc */
11829566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &submats[0]));
11839566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(submats[0], 0, 0, PETSC_DETERMINE, PETSC_DETERMINE));
11849566063dSJacob Faibussowitsch       PetscCall(MatSetType(submats[0], MATDUMMY));
11854698041cSHong Zhang 
11865c39f6d9SHong Zhang       /* create struct Mat_SubSppt and attached it to submat */
11874dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&smat_i));
11884698041cSHong Zhang       submats[0]->data = (void *)smat_i;
11894698041cSHong Zhang 
11904698041cSHong Zhang       smat_i->destroy          = submats[0]->ops->destroy;
1191f68bb481SHong Zhang       submats[0]->ops->destroy = MatDestroySubMatrix_Dummy;
11924698041cSHong Zhang       submats[0]->factortype   = C->factortype;
11934698041cSHong Zhang 
1194006cef31SHong Zhang       smat_i->id          = 0;
11954698041cSHong Zhang       smat_i->nrqs        = nrqs;
11964698041cSHong Zhang       smat_i->nrqr        = nrqr;
11974698041cSHong Zhang       smat_i->rbuf1       = rbuf1;
11984698041cSHong Zhang       smat_i->rbuf2       = rbuf2;
11994698041cSHong Zhang       smat_i->rbuf3       = rbuf3;
12004698041cSHong Zhang       smat_i->sbuf2       = sbuf2;
12014698041cSHong Zhang       smat_i->req_source2 = req_source2;
12024698041cSHong Zhang 
12034698041cSHong Zhang       smat_i->sbuf1 = sbuf1;
12044698041cSHong Zhang       smat_i->ptr   = ptr;
12054698041cSHong Zhang       smat_i->tmp   = tmp;
12064698041cSHong Zhang       smat_i->ctr   = ctr;
12074698041cSHong Zhang 
12084698041cSHong Zhang       smat_i->pa          = pa;
12094698041cSHong Zhang       smat_i->req_size    = req_size;
12104698041cSHong Zhang       smat_i->req_source1 = req_source1;
12114698041cSHong Zhang 
1212c67fe082SHong Zhang       smat_i->allcolumns = PETSC_FALSE;
12134698041cSHong Zhang       smat_i->singleis   = PETSC_FALSE;
12144698041cSHong Zhang       smat_i->row2proc   = NULL;
12154698041cSHong Zhang       smat_i->rmap       = NULL;
12164698041cSHong Zhang       smat_i->cmap       = NULL;
12174698041cSHong Zhang     }
12184698041cSHong Zhang 
12199566063dSJacob Faibussowitsch     if (ismax) PetscCall(PetscFree(lens[0]));
12209566063dSJacob Faibussowitsch     PetscCall(PetscFree(lens));
1221175c2bc5SJunchao Zhang     if (sbuf_aj) {
12229566063dSJacob Faibussowitsch       PetscCall(PetscFree(sbuf_aj[0]));
12239566063dSJacob Faibussowitsch       PetscCall(PetscFree(sbuf_aj));
1224175c2bc5SJunchao Zhang     }
122517df9f7cSHong Zhang 
122617df9f7cSHong Zhang   } /* endof scall == MAT_INITIAL_MATRIX */
122717df9f7cSHong Zhang 
122817df9f7cSHong Zhang   /* Post recv matrix values */
1229bbc89d27SHong Zhang   if (!ijonly) {
12309566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag4));
12319566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &rbuf4));
12329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqs, &r_waits4));
123317df9f7cSHong Zhang     for (i = 0; i < nrqs; ++i) {
12349566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(rbuf2[i][0] * bs2, &rbuf4[i]));
12359566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(rbuf4[i], rbuf2[i][0] * bs2, MPIU_SCALAR, req_source2[i], tag4, comm, r_waits4 + i));
123617df9f7cSHong Zhang     }
123717df9f7cSHong Zhang 
123817df9f7cSHong Zhang     /* Allocate sending buffers for a->a, and send them off */
12399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &sbuf_aa));
124017df9f7cSHong Zhang     for (i = 0, j = 0; i < nrqr; i++) j += req_size[i];
12419e686edcSHong Zhang 
12429566063dSJacob Faibussowitsch     if (nrqr) PetscCall(PetscMalloc1(j * bs2, &sbuf_aa[0]));
1243a42ab0d6SHong Zhang     for (i = 1; i < nrqr; i++) sbuf_aa[i] = sbuf_aa[i - 1] + req_size[i - 1] * bs2;
124417df9f7cSHong Zhang 
12459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrqr, &s_waits4));
124617df9f7cSHong Zhang 
124717df9f7cSHong Zhang     for (i = 0; i < nrqr; i++) {
124817df9f7cSHong Zhang       rbuf1_i   = rbuf1[i];
124917df9f7cSHong Zhang       sbuf_aa_i = sbuf_aa[i];
125017df9f7cSHong Zhang       ct1       = 2 * rbuf1_i[0] + 1;
125117df9f7cSHong Zhang       ct2       = 0;
125217df9f7cSHong Zhang       for (j = 1, max1 = rbuf1_i[0]; j <= max1; j++) {
125317df9f7cSHong Zhang         kmax = rbuf1_i[2 * j];
125417df9f7cSHong Zhang         for (k = 0; k < kmax; k++, ct1++) {
125517df9f7cSHong Zhang           row    = rbuf1_i[ct1] - rstart;
1256a42ab0d6SHong Zhang           nzA    = a_i[row + 1] - a_i[row];
1257a42ab0d6SHong Zhang           nzB    = b_i[row + 1] - b_i[row];
125817df9f7cSHong Zhang           ncols  = nzA + nzB;
12598e3a54c0SPierre Jolivet           cworkB = PetscSafePointerPlusOffset(b_j, b_i[row]);
12608e3a54c0SPierre Jolivet           vworkA = PetscSafePointerPlusOffset(a_a, a_i[row] * bs2);
12618e3a54c0SPierre Jolivet           vworkB = PetscSafePointerPlusOffset(b_a, b_i[row] * bs2);
126217df9f7cSHong Zhang 
126317df9f7cSHong Zhang           /* load the column values for this row into vals*/
1264a42ab0d6SHong Zhang           vals = sbuf_aa_i + ct2 * bs2;
1265a42ab0d6SHong Zhang           for (l = 0; l < nzB; l++) {
1266a42ab0d6SHong Zhang             if ((bmap[cworkB[l]]) < cstart) {
12679566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(vals + l * bs2, vworkB + l * bs2, bs2));
1268a42ab0d6SHong Zhang             } else break;
1269a42ab0d6SHong Zhang           }
1270a42ab0d6SHong Zhang           imark = l;
127148a46eb9SPierre Jolivet           for (l = 0; l < nzA; l++) PetscCall(PetscArraycpy(vals + (imark + l) * bs2, vworkA + l * bs2, bs2));
127248a46eb9SPierre Jolivet           for (l = imark; l < nzB; l++) PetscCall(PetscArraycpy(vals + (nzA + l) * bs2, vworkB + l * bs2, bs2));
127317df9f7cSHong Zhang 
127417df9f7cSHong Zhang           ct2 += ncols;
127517df9f7cSHong Zhang         }
127617df9f7cSHong Zhang       }
12779566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sbuf_aa_i, req_size[i] * bs2, MPIU_SCALAR, req_source1[i], tag4, comm, s_waits4 + i));
127817df9f7cSHong Zhang     }
1279bbc89d27SHong Zhang   }
128017df9f7cSHong Zhang 
128117df9f7cSHong Zhang   /* Assemble the matrices */
128217df9f7cSHong Zhang   /* First assemble the local rows */
128317df9f7cSHong Zhang   for (i = 0; i < ismax; i++) {
128417df9f7cSHong Zhang     row2proc_i = row2proc[i];
128517df9f7cSHong Zhang     subc       = (Mat_SeqBAIJ *)submats[i]->data;
128617df9f7cSHong Zhang     imat_ilen  = subc->ilen;
128717df9f7cSHong Zhang     imat_j     = subc->j;
128817df9f7cSHong Zhang     imat_i     = subc->i;
128917df9f7cSHong Zhang     imat_a     = subc->a;
129017df9f7cSHong Zhang 
129117df9f7cSHong Zhang     if (!allcolumns[i]) cmap_i = cmap[i];
129217df9f7cSHong Zhang     rmap_i = rmap[i];
129317df9f7cSHong Zhang     irow_i = irow[i];
129417df9f7cSHong Zhang     jmax   = nrow[i];
129517df9f7cSHong Zhang     for (j = 0; j < jmax; j++) {
1296ec3c739cSHong Zhang       if (allrows[i]) row = j;
1297ec3c739cSHong Zhang       else row = irow_i[j];
129817df9f7cSHong Zhang       proc = row2proc_i[j];
1299a42ab0d6SHong Zhang 
130017df9f7cSHong Zhang       if (proc == rank) {
1301a42ab0d6SHong Zhang         row    = row - rstart;
1302a42ab0d6SHong Zhang         nzA    = a_i[row + 1] - a_i[row];
1303a42ab0d6SHong Zhang         nzB    = b_i[row + 1] - b_i[row];
1304a42ab0d6SHong Zhang         cworkA = a_j + a_i[row];
13058e3a54c0SPierre Jolivet         cworkB = PetscSafePointerPlusOffset(b_j, b_i[row]);
1306a42ab0d6SHong Zhang         if (!ijonly) {
1307a42ab0d6SHong Zhang           vworkA = a_a + a_i[row] * bs2;
13088e3a54c0SPierre Jolivet           vworkB = PetscSafePointerPlusOffset(b_a, b_i[row] * bs2);
1309a42ab0d6SHong Zhang         }
1310059f6b74SHong Zhang 
1311059f6b74SHong Zhang         if (allrows[i]) {
1312059f6b74SHong Zhang           row = row + rstart;
1313059f6b74SHong Zhang         } else {
1314a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1315eec179cfSJacob Faibussowitsch           PetscCall(PetscHMapIGetWithDefault(rmap_i, row + rstart + 1, 0, &row));
1316a42ab0d6SHong Zhang           row--;
1317121971b2SHong Zhang 
131808401ef6SPierre Jolivet           PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table");
1319a42ab0d6SHong Zhang #else
1320a42ab0d6SHong Zhang           row = rmap_i[row + rstart];
1321a42ab0d6SHong Zhang #endif
1322059f6b74SHong Zhang         }
1323a42ab0d6SHong Zhang         mat_i = imat_i[row];
13248e3a54c0SPierre Jolivet         if (!ijonly) mat_a = PetscSafePointerPlusOffset(imat_a, mat_i * bs2);
13258e3a54c0SPierre Jolivet         mat_j = PetscSafePointerPlusOffset(imat_j, mat_i);
132680d31651SHong Zhang         ilen  = imat_ilen[row];
1327a42ab0d6SHong Zhang 
1328a42ab0d6SHong Zhang         /* load the column indices for this row into cols*/
1329a42ab0d6SHong Zhang         if (!allcolumns[i]) {
1330a42ab0d6SHong Zhang           for (l = 0; l < nzB; l++) {
1331a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1332a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1333eec179cfSJacob Faibussowitsch               PetscCall(PetscHMapIGetWithDefault(cmap_i, ctmp + 1, 0, &tcol));
1334a42ab0d6SHong Zhang               if (tcol) {
1335a42ab0d6SHong Zhang #else
1336a42ab0d6SHong Zhang               if ((tcol = cmap_i[ctmp])) {
1337a42ab0d6SHong Zhang #endif
1338a42ab0d6SHong Zhang                 *mat_j++ = tcol - 1;
13399566063dSJacob Faibussowitsch                 PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2));
1340a42ab0d6SHong Zhang                 mat_a += bs2;
134180d31651SHong Zhang                 ilen++;
1342a42ab0d6SHong Zhang               }
1343a42ab0d6SHong Zhang             } else break;
1344a42ab0d6SHong Zhang           }
1345a42ab0d6SHong Zhang           imark = l;
1346a42ab0d6SHong Zhang           for (l = 0; l < nzA; l++) {
1347a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1348eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapIGetWithDefault(cmap_i, cstart + cworkA[l] + 1, 0, &tcol));
1349a42ab0d6SHong Zhang             if (tcol) {
1350a42ab0d6SHong Zhang #else
1351a42ab0d6SHong Zhang             if ((tcol = cmap_i[cstart + cworkA[l]])) {
1352a42ab0d6SHong Zhang #endif
1353a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1354a42ab0d6SHong Zhang               if (!ijonly) {
13559566063dSJacob Faibussowitsch                 PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2));
1356a42ab0d6SHong Zhang                 mat_a += bs2;
1357a42ab0d6SHong Zhang               }
135880d31651SHong Zhang               ilen++;
1359a42ab0d6SHong Zhang             }
1360a42ab0d6SHong Zhang           }
1361a42ab0d6SHong Zhang           for (l = imark; l < nzB; l++) {
1362a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1363eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapIGetWithDefault(cmap_i, bmap[cworkB[l]] + 1, 0, &tcol));
1364a42ab0d6SHong Zhang             if (tcol) {
1365a42ab0d6SHong Zhang #else
1366a42ab0d6SHong Zhang             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1367a42ab0d6SHong Zhang #endif
1368a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1369a42ab0d6SHong Zhang               if (!ijonly) {
13709566063dSJacob Faibussowitsch                 PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2));
1371a42ab0d6SHong Zhang                 mat_a += bs2;
1372a42ab0d6SHong Zhang               }
137380d31651SHong Zhang               ilen++;
1374a42ab0d6SHong Zhang             }
1375a42ab0d6SHong Zhang           }
1376a42ab0d6SHong Zhang         } else { /* allcolumns */
1377a42ab0d6SHong Zhang           for (l = 0; l < nzB; l++) {
1378a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1379a42ab0d6SHong Zhang               *mat_j++ = ctmp;
13809566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2));
1381a42ab0d6SHong Zhang               mat_a += bs2;
138280d31651SHong Zhang               ilen++;
1383a42ab0d6SHong Zhang             } else break;
1384a42ab0d6SHong Zhang           }
1385a42ab0d6SHong Zhang           imark = l;
1386a42ab0d6SHong Zhang           for (l = 0; l < nzA; l++) {
1387a42ab0d6SHong Zhang             *mat_j++ = cstart + cworkA[l];
1388a42ab0d6SHong Zhang             if (!ijonly) {
13899566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2));
1390a42ab0d6SHong Zhang               mat_a += bs2;
1391a42ab0d6SHong Zhang             }
139280d31651SHong Zhang             ilen++;
1393a42ab0d6SHong Zhang           }
1394a42ab0d6SHong Zhang           for (l = imark; l < nzB; l++) {
1395a42ab0d6SHong Zhang             *mat_j++ = bmap[cworkB[l]];
1396a42ab0d6SHong Zhang             if (!ijonly) {
13979566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2));
1398a42ab0d6SHong Zhang               mat_a += bs2;
1399a42ab0d6SHong Zhang             }
140080d31651SHong Zhang             ilen++;
1401a42ab0d6SHong Zhang           }
1402a42ab0d6SHong Zhang         }
140380d31651SHong Zhang         imat_ilen[row] = ilen;
140417df9f7cSHong Zhang       }
140517df9f7cSHong Zhang     }
140617df9f7cSHong Zhang   }
140717df9f7cSHong Zhang 
140817df9f7cSHong Zhang   /* Now assemble the off proc rows */
140948a46eb9SPierre Jolivet   if (!ijonly) PetscCallMPI(MPI_Waitall(nrqs, r_waits4, MPI_STATUSES_IGNORE));
141017df9f7cSHong Zhang   for (tmp2 = 0; tmp2 < nrqs; tmp2++) {
141117df9f7cSHong Zhang     sbuf1_i = sbuf1[pa[tmp2]];
141217df9f7cSHong Zhang     jmax    = sbuf1_i[0];
141317df9f7cSHong Zhang     ct1     = 2 * jmax + 1;
141417df9f7cSHong Zhang     ct2     = 0;
141517df9f7cSHong Zhang     rbuf2_i = rbuf2[tmp2];
141617df9f7cSHong Zhang     rbuf3_i = rbuf3[tmp2];
14175dd0c08cSHong Zhang     if (!ijonly) rbuf4_i = rbuf4[tmp2];
141817df9f7cSHong Zhang     for (j = 1; j <= jmax; j++) {
141917df9f7cSHong Zhang       is_no  = sbuf1_i[2 * j - 1];
142017df9f7cSHong Zhang       rmap_i = rmap[is_no];
142117df9f7cSHong Zhang       if (!allcolumns[is_no]) cmap_i = cmap[is_no];
142217df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ *)submats[is_no]->data;
142317df9f7cSHong Zhang       imat_ilen = subc->ilen;
142417df9f7cSHong Zhang       imat_j    = subc->j;
142517df9f7cSHong Zhang       imat_i    = subc->i;
14265dd0c08cSHong Zhang       if (!ijonly) imat_a = subc->a;
142717df9f7cSHong Zhang       max1 = sbuf1_i[2 * j];
14285dd0c08cSHong Zhang       for (k = 0; k < max1; k++, ct1++) { /* for each recved block row */
142917df9f7cSHong Zhang         row = sbuf1_i[ct1];
1430059f6b74SHong Zhang 
1431059f6b74SHong Zhang         if (allrows[is_no]) {
1432059f6b74SHong Zhang           row = sbuf1_i[ct1];
1433059f6b74SHong Zhang         } else {
143417df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1435eec179cfSJacob Faibussowitsch           PetscCall(PetscHMapIGetWithDefault(rmap_i, row + 1, 0, &row));
143617df9f7cSHong Zhang           row--;
143708401ef6SPierre Jolivet           PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table");
143817df9f7cSHong Zhang #else
143917df9f7cSHong Zhang           row = rmap_i[row];
144017df9f7cSHong Zhang #endif
1441059f6b74SHong Zhang         }
144217df9f7cSHong Zhang         ilen  = imat_ilen[row];
144317df9f7cSHong Zhang         mat_i = imat_i[row];
14445dd0c08cSHong Zhang         if (!ijonly) mat_a = imat_a + mat_i * bs2;
144517df9f7cSHong Zhang         mat_j = imat_j + mat_i;
144617df9f7cSHong Zhang         max2  = rbuf2_i[ct1];
144717df9f7cSHong Zhang         if (!allcolumns[is_no]) {
144817df9f7cSHong Zhang           for (l = 0; l < max2; l++, ct2++) {
144917df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1450eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapIGetWithDefault(cmap_i, rbuf3_i[ct2] + 1, 0, &tcol));
145117df9f7cSHong Zhang #else
145217df9f7cSHong Zhang             tcol = cmap_i[rbuf3_i[ct2]];
145317df9f7cSHong Zhang #endif
145417df9f7cSHong Zhang             if (tcol) {
145517df9f7cSHong Zhang               *mat_j++ = tcol - 1;
14565dd0c08cSHong Zhang               if (!ijonly) {
14579566063dSJacob Faibussowitsch                 PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2));
14585dd0c08cSHong Zhang                 mat_a += bs2;
14595dd0c08cSHong Zhang               }
146017df9f7cSHong Zhang               ilen++;
146117df9f7cSHong Zhang             }
146217df9f7cSHong Zhang           }
146317df9f7cSHong Zhang         } else { /* allcolumns */
146417df9f7cSHong Zhang           for (l = 0; l < max2; l++, ct2++) {
146517df9f7cSHong Zhang             *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */
14665dd0c08cSHong Zhang             if (!ijonly) {
14679566063dSJacob Faibussowitsch               PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2));
14685dd0c08cSHong Zhang               mat_a += bs2;
14695dd0c08cSHong Zhang             }
147017df9f7cSHong Zhang             ilen++;
147117df9f7cSHong Zhang           }
147217df9f7cSHong Zhang         }
147317df9f7cSHong Zhang         imat_ilen[row] = ilen;
147417df9f7cSHong Zhang       }
147517df9f7cSHong Zhang     }
147617df9f7cSHong Zhang   }
147717df9f7cSHong Zhang 
147880d31651SHong Zhang   if (!iscsorted) { /* sort column indices of the rows */
147980d31651SHong Zhang     MatScalar *work;
148080d31651SHong Zhang 
14819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs2, &work));
148217df9f7cSHong Zhang     for (i = 0; i < ismax; i++) {
148317df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ *)submats[i]->data;
148480d31651SHong Zhang       imat_ilen = subc->ilen;
148517df9f7cSHong Zhang       imat_j    = subc->j;
148617df9f7cSHong Zhang       imat_i    = subc->i;
14874b6ceb0dSHong Zhang       if (!ijonly) imat_a = subc->a;
148817df9f7cSHong Zhang       if (allcolumns[i]) continue;
14894b6ceb0dSHong Zhang 
149017df9f7cSHong Zhang       jmax = nrow[i];
149117df9f7cSHong Zhang       for (j = 0; j < jmax; j++) {
149280d31651SHong Zhang         mat_i = imat_i[j];
149380d31651SHong Zhang         mat_j = imat_j + mat_i;
14944b6ceb0dSHong Zhang         ilen  = imat_ilen[j];
149580d31651SHong Zhang         if (ijonly) {
14969566063dSJacob Faibussowitsch           PetscCall(PetscSortInt(ilen, mat_j));
149780d31651SHong Zhang         } else {
14984b6ceb0dSHong Zhang           mat_a = imat_a + mat_i * bs2;
14999566063dSJacob Faibussowitsch           PetscCall(PetscSortIntWithDataArray(ilen, mat_j, mat_a, bs2 * sizeof(MatScalar), work));
150017df9f7cSHong Zhang         }
150117df9f7cSHong Zhang       }
150217df9f7cSHong Zhang     }
15039566063dSJacob Faibussowitsch     PetscCall(PetscFree(work));
150480d31651SHong Zhang   }
150517df9f7cSHong Zhang 
1506bbc89d27SHong Zhang   if (!ijonly) {
15079566063dSJacob Faibussowitsch     PetscCall(PetscFree(r_waits4));
15089566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrqr, s_waits4, MPI_STATUSES_IGNORE));
15099566063dSJacob Faibussowitsch     PetscCall(PetscFree(s_waits4));
1510bbc89d27SHong Zhang   }
151117df9f7cSHong Zhang 
151217df9f7cSHong Zhang   /* Restore the indices */
151317df9f7cSHong Zhang   for (i = 0; i < ismax; i++) {
151448a46eb9SPierre Jolivet     if (!allrows[i]) PetscCall(ISRestoreIndices(isrow[i], irow + i));
151548a46eb9SPierre Jolivet     if (!allcolumns[i]) PetscCall(ISRestoreIndices(iscol[i], icol + i));
151617df9f7cSHong Zhang   }
151717df9f7cSHong Zhang 
151817df9f7cSHong Zhang   for (i = 0; i < ismax; i++) {
15199566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY));
15209566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY));
152117df9f7cSHong Zhang   }
152217df9f7cSHong Zhang 
15239566063dSJacob Faibussowitsch   PetscCall(PetscFree5(*(PetscInt ***)&irow, *(PetscInt ***)&icol, nrow, ncol, issorted));
15249566063dSJacob Faibussowitsch   PetscCall(PetscFree5(row2proc, cmap, rmap, allcolumns, allrows));
152517df9f7cSHong Zhang 
1526bbc89d27SHong Zhang   if (!ijonly) {
1527175c2bc5SJunchao Zhang     if (sbuf_aa) {
15289566063dSJacob Faibussowitsch       PetscCall(PetscFree(sbuf_aa[0]));
15299566063dSJacob Faibussowitsch       PetscCall(PetscFree(sbuf_aa));
1530175c2bc5SJunchao Zhang     }
153117df9f7cSHong Zhang 
153248a46eb9SPierre Jolivet     for (i = 0; i < nrqs; ++i) PetscCall(PetscFree(rbuf4[i]));
15339566063dSJacob Faibussowitsch     PetscCall(PetscFree(rbuf4));
1534bbc89d27SHong Zhang   }
15357a868f3eSHong Zhang   c->ijonly = PETSC_FALSE; /* set back to the default */
15363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1537d5b485abSSatish Balay }
1538