xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
14b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov)
15d5b485abSSatish Balay {
16d0f46423SBarry Smith   PetscInt       i,N=C->cmap->N, bs=C->rmap->bs;
1736f4e84dSSatish Balay   IS             *is_new;
1836f4e84dSSatish Balay 
193a40ed3dSBarry Smith   PetscFunctionBegin;
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(imax,&is_new));
21df36731bSBarry Smith   /* Convert the indices into block format */
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new));
232c71b3e2SJacob Faibussowitsch   PetscCheckFalse(ov < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
24d5b485abSSatish Balay   for (i=0; i<ov; ++i) {
25*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new));
26d5b485abSSatish Balay   }
27*5f80ce2aSJacob Faibussowitsch   for (i=0; i<imax; i++) CHKERRQ(ISDestroy(&is[i]));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISExpandIndicesGeneral(N,N,bs,imax,is_new,is));
29*5f80ce2aSJacob Faibussowitsch   for (i=0; i<imax; i++) CHKERRQ(ISDestroy(&is_new[i]));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is_new));
313a40ed3dSBarry Smith   PetscFunctionReturn(0);
32d5b485abSSatish Balay }
33d5b485abSSatish Balay 
34d5b485abSSatish Balay /*
35d5b485abSSatish Balay   Sample message format:
36d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
370e9b0e7eSHong Zhang   to index sets is[1], is[5]
38d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
39d5b485abSSatish Balay   -----------
40d5b485abSSatish Balay   mesg [1] = 1 => is[1]
41d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
42d5b485abSSatish Balay   -----------
43d5b485abSSatish Balay   mesg [5] = 5  => is[5]
44d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
45d5b485abSSatish Balay   -----------
46d5b485abSSatish Balay   mesg [7]
470e9b0e7eSHong Zhang   mesg [n]  data(is[1])
48d5b485abSSatish Balay   -----------
49d5b485abSSatish Balay   mesg[n+1]
50d5b485abSSatish Balay   mesg[m]  data(is[5])
51d5b485abSSatish Balay   -----------
52d5b485abSSatish Balay 
53d5b485abSSatish Balay   Notes:
54d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
552d4ee042Sprj-   nrqr - no of requests received (which have to be or which have been processed)
56d5b485abSSatish Balay */
57db41cccfSHong Zhang PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[])
58d5b485abSSatish Balay {
59df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
605d0c19d7SBarry Smith   const PetscInt **idx,*idx_i;
6124c049a4SHong Zhang   PetscInt       *n,*w3,*w4,**data,len;
62b24ad042SBarry Smith   PetscMPIInt    size,rank,tag1,tag2,*w2,*w1,nrqr;
63131c27b5Sprj-   PetscInt       Mbs,i,j,k,**rbuf,row,nrqs,msz,**outdat,**ptr;
648f9f447aSHong Zhang   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p;
65131c27b5Sprj-   PetscMPIInt    *onodes1,*olengths1,*onodes2,*olengths2,proc=-1;
66f1af5d2fSBarry Smith   PetscBT        *table;
67749130a8SStefano Zampini   MPI_Comm       comm,*iscomms;
68d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
698f9f447aSHong Zhang   char           *t_p;
70d5b485abSSatish Balay 
713a40ed3dSBarry Smith   PetscFunctionBegin;
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)C,&comm));
73d5b485abSSatish Balay   size = c->size;
74d5b485abSSatish Balay   rank = c->rank;
75df36731bSBarry Smith   Mbs  = c->Mbs;
76d5b485abSSatish Balay 
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetNewTag((PetscObject)C,&tag1));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetNewTag((PetscObject)C,&tag2));
79c7dd2924SSatish Balay 
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(imax+1,(PetscInt***)&idx,imax,&n));
8124c049a4SHong Zhang 
82d5b485abSSatish Balay   for (i=0; i<imax; i++) {
83*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(is[i],&idx[i]));
84*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(is[i],&n[i]));
85d5b485abSSatish Balay   }
86d5b485abSSatish Balay 
87d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
88d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4));
90d5b485abSSatish Balay   for (i=0; i<imax; i++) {
91*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(w4,size)); /* initialise work vector*/
92d5b485abSSatish Balay     idx_i = idx[i];
93d5b485abSSatish Balay     len   = n[i];
94d5b485abSSatish Balay     for (j=0; j<len; j++) {
95d5b485abSSatish Balay       row = idx_i[j];
962c71b3e2SJacob Faibussowitsch       PetscCheckFalse(row < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries");
97*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc));
98d5b485abSSatish Balay       w4[proc]++;
99d5b485abSSatish Balay     }
100d5b485abSSatish Balay     for (j=0; j<size; j++) {
101d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
102d5b485abSSatish Balay     }
103d5b485abSSatish Balay   }
104d5b485abSSatish Balay 
105d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
106d5b485abSSatish Balay   msz      = 0;              /* total mesg length (for all proc */
1070e9b0e7eSHong Zhang   w1[rank] = 0;              /* no mesg sent to itself */
108d5b485abSSatish Balay   w3[rank] = 0;
109d5b485abSSatish Balay   for (i=0; i<size; i++) {
110d5b485abSSatish Balay     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
111d5b485abSSatish Balay   }
112d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrqs,&pa));
114d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
115d5b485abSSatish Balay     if (w1[i]) {pa[j] = i; j++;}
116d5b485abSSatish Balay   }
117d5b485abSSatish Balay 
118d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
119d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
120d5b485abSSatish Balay     j      = pa[i];
121d5b485abSSatish Balay     w1[j] += w2[j] + 2*w3[j];
122d5b485abSSatish Balay     msz   += w1[j];
123d5b485abSSatish Balay   }
124d5b485abSSatish Balay 
125c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscGatherNumberOfMessages(comm,w2,w1,&nrqr));
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1));
128d5b485abSSatish Balay 
129c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1));
131d5b485abSSatish Balay 
132d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr));
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(outdat,size));
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(ptr,size));
136d5b485abSSatish Balay   {
137b24ad042SBarry Smith     PetscInt *iptr = tmp,ict  = 0;
138d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
139d5b485abSSatish Balay       j         = pa[i];
140d5b485abSSatish Balay       iptr     +=  ict;
141d5b485abSSatish Balay       outdat[j] = iptr;
142d5b485abSSatish Balay       ict       = w1[j];
143d5b485abSSatish Balay     }
144d5b485abSSatish Balay   }
145d5b485abSSatish Balay 
146d5b485abSSatish Balay   /* Form the outgoing messages */
147d5b485abSSatish Balay   /*plug in the headers*/
148d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
149d5b485abSSatish Balay     j            = pa[i];
150d5b485abSSatish Balay     outdat[j][0] = 0;
151*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(outdat[j]+1,2*w3[j]));
152d5b485abSSatish Balay     ptr[j]       = outdat[j] + 2*w3[j] + 1;
153d5b485abSSatish Balay   }
154d5b485abSSatish Balay 
155d5b485abSSatish Balay   /* Memory for doing local proc's work*/
156d5b485abSSatish Balay   {
157*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mbs*imax,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,&t_p));
158d5b485abSSatish Balay 
159bbd702dbSSatish Balay     for (i=0; i<imax; i++) {
160b6410449SSatish Balay       table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
161bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
162d5b485abSSatish Balay     }
163d5b485abSSatish Balay   }
164d5b485abSSatish Balay 
165d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
166d5b485abSSatish Balay   {
167b24ad042SBarry Smith     PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j;
168f1af5d2fSBarry Smith     PetscBT  table_i;
169d5b485abSSatish Balay 
170d5b485abSSatish Balay     for (i=0; i<imax; i++) {
171*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctr,size));
172d5b485abSSatish Balay       n_i     = n[i];
173d5b485abSSatish Balay       table_i = table[i];
174d5b485abSSatish Balay       idx_i   = idx[i];
175d5b485abSSatish Balay       data_i  = data[i];
176d5b485abSSatish Balay       isz_i   = isz[i];
177d5b485abSSatish Balay       for (j=0; j<n_i; j++) {   /* parse the indices of each IS */
178d5b485abSSatish Balay         row  = idx_i[j];
179*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc));
180d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
181d5b485abSSatish Balay           ctr[proc]++;
182d5b485abSSatish Balay           *ptr[proc] = row;
183d5b485abSSatish Balay           ptr[proc]++;
184d6b45a43SBarry Smith         } else { /* Update the local table */
18526fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row;
186d5b485abSSatish Balay         }
187d5b485abSSatish Balay       }
188d5b485abSSatish Balay       /* Update the headers for the current IS */
189d5b485abSSatish Balay       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
190d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
191d5b485abSSatish Balay           outdat_j        = outdat[j];
192d5b485abSSatish Balay           k               = ++outdat_j[0];
193d5b485abSSatish Balay           outdat_j[2*k]   = ctr_j;
194d5b485abSSatish Balay           outdat_j[2*k-1] = i;
195d5b485abSSatish Balay         }
196d5b485abSSatish Balay       }
197d5b485abSSatish Balay       isz[i] = isz_i;
198d5b485abSSatish Balay     }
199d5b485abSSatish Balay   }
200d5b485abSSatish Balay 
201d5b485abSSatish Balay   /*  Now  post the sends */
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrqs,&s_waits1));
203d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
204d5b485abSSatish Balay     j    = pa[i];
205*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i));
206d5b485abSSatish Balay   }
207d5b485abSSatish Balay 
208d5b485abSSatish Balay   /* No longer need the original indices*/
209d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
210*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(is[i],idx+i));
211d5b485abSSatish Balay   }
212*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(*(PetscInt***)&idx,n));
213d5b485abSSatish Balay 
214*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(imax,&iscomms));
215d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
216*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]),&iscomms[i],NULL));
217*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is[i]));
218d5b485abSSatish Balay   }
219d5b485abSSatish Balay 
220d5b485abSSatish Balay   /* Do Local work*/
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data));
222d5b485abSSatish Balay 
223d5b485abSSatish Balay   /* Receive messages*/
224*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(nrqr,r_waits1,MPI_STATUSES_IGNORE));
225*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(nrqs,s_waits1,MPI_STATUSES_IGNORE));
226d5b485abSSatish Balay 
227d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
228*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(outdat,ptr,tmp,ctr));
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(w1,w2,w3,w4));
230d5b485abSSatish Balay 
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrqr,&xdata));
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrqr,&isz1));
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1));
234175c2bc5SJunchao Zhang   if (rbuf) {
235*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rbuf[0]));
236*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rbuf));
237175c2bc5SJunchao Zhang   }
238d5b485abSSatish Balay 
239d5b485abSSatish Balay   /* Send the data back*/
240d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
241d5b485abSSatish Balay   {
242b24ad042SBarry Smith     PetscMPIInt *rw1;
243d5b485abSSatish Balay 
244*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(size,&rw1));
245c7dd2924SSatish Balay 
246d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
247175c2bc5SJunchao Zhang       proc = onodes1[i];
248d5b485abSSatish Balay       rw1[proc] = isz1[i];
249d5b485abSSatish Balay     }
250d5b485abSSatish Balay 
251c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
252*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2));
253*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rw1));
254c7dd2924SSatish Balay   }
255c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
256*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2));
257d5b485abSSatish Balay 
258d5b485abSSatish Balay   /*  Now  post the sends */
259*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrqr,&s_waits2));
260d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
261175c2bc5SJunchao Zhang     j    = onodes1[i];
262*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i));
263d5b485abSSatish Balay   }
264d5b485abSSatish Balay 
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(onodes1));
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(olengths1));
267175c2bc5SJunchao Zhang 
268d5b485abSSatish Balay   /* receive work done on other processors*/
269d5b485abSSatish Balay   {
270b24ad042SBarry Smith     PetscMPIInt idex;
271b24ad042SBarry Smith     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
272f1af5d2fSBarry Smith     PetscBT     table_i;
273d5b485abSSatish Balay 
274d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
275*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Waitany(nrqs,r_waits2,&idex,MPI_STATUS_IGNORE));
276d5b485abSSatish Balay       /* Process the message*/
277999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
278d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
279999d9058SBarry Smith       jmax    = rbuf2[idex][0];
280d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
281d5b485abSSatish Balay         max     = rbuf2_i[2*j];
282d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
283d5b485abSSatish Balay         isz_i   = isz[is_no];
284d5b485abSSatish Balay         data_i  = data[is_no];
285d5b485abSSatish Balay         table_i = table[is_no];
286d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
287d5b485abSSatish Balay           row = rbuf2_i[ct1];
28826fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row;
289d5b485abSSatish Balay         }
290d5b485abSSatish Balay         isz[is_no] = isz_i;
291d5b485abSSatish Balay       }
292d5b485abSSatish Balay     }
293*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrqr,s_waits2,MPI_STATUSES_IGNORE));
294d5b485abSSatish Balay   }
295d5b485abSSatish Balay 
296d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
297*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(iscomms[i],isz[i],data[i],PETSC_COPY_VALUES,is+i));
298*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCommDestroy(&iscomms[i]));
299d5b485abSSatish Balay   }
300d5b485abSSatish Balay 
301*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(iscomms));
302*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(onodes2));
303*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(olengths2));
304c7dd2924SSatish Balay 
305*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pa));
306175c2bc5SJunchao Zhang   if (rbuf2) {
307*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rbuf2[0]));
308*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rbuf2));
309175c2bc5SJunchao Zhang   }
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(s_waits1));
311*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(r_waits1));
312*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(s_waits2));
313*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(r_waits2));
314*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree5(table,data,isz,d_p,t_p));
315175c2bc5SJunchao Zhang   if (xdata) {
316*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(xdata[0]));
317*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(xdata));
318175c2bc5SJunchao Zhang   }
319*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(isz1));
3203a40ed3dSBarry Smith   PetscFunctionReturn(0);
321d5b485abSSatish Balay }
322d5b485abSSatish Balay 
323d5b485abSSatish Balay /*
324df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
325d5b485abSSatish Balay        the work on the local processor.
326d5b485abSSatish Balay 
327d5b485abSSatish Balay      Inputs:
328df36731bSBarry Smith       C      - MAT_MPIBAIJ;
329d5b485abSSatish Balay       imax - total no of index sets processed at a time;
330df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
331d5b485abSSatish Balay 
332d5b485abSSatish Balay      Output:
3330e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
334d5b485abSSatish Balay                to each index set;
335d5b485abSSatish Balay       data   - pointer to the solutions
336d5b485abSSatish Balay */
337b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data)
338d5b485abSSatish Balay {
339df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
340d5b485abSSatish Balay   Mat         A  = c->A,B = c->B;
341df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
342b24ad042SBarry Smith   PetscInt    start,end,val,max,rstart,cstart,*ai,*aj;
343b24ad042SBarry Smith   PetscInt    *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
344f1af5d2fSBarry Smith   PetscBT     table_i;
345d5b485abSSatish Balay 
3463a40ed3dSBarry Smith   PetscFunctionBegin;
347899cda47SBarry Smith   rstart = c->rstartbs;
348899cda47SBarry Smith   cstart = c->cstartbs;
349d5b485abSSatish Balay   ai     = a->i;
350df36731bSBarry Smith   aj     = a->j;
351d5b485abSSatish Balay   bi     = b->i;
352df36731bSBarry Smith   bj     = b->j;
353d5b485abSSatish Balay   garray = c->garray;
354d5b485abSSatish Balay 
355d5b485abSSatish Balay   for (i=0; i<imax; i++) {
356d5b485abSSatish Balay     data_i  = data[i];
357d5b485abSSatish Balay     table_i = table[i];
358d5b485abSSatish Balay     isz_i   = isz[i];
359d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
360d5b485abSSatish Balay       row   = data_i[j] - rstart;
361d5b485abSSatish Balay       start = ai[row];
362d5b485abSSatish Balay       end   = ai[row+1];
363d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
364df36731bSBarry Smith         val = aj[k] + cstart;
36526fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val;
366d5b485abSSatish Balay       }
367d5b485abSSatish Balay       start = bi[row];
368d5b485abSSatish Balay       end   = bi[row+1];
369d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
370df36731bSBarry Smith         val = garray[bj[k]];
37126fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val;
372d5b485abSSatish Balay       }
373d5b485abSSatish Balay     }
374d5b485abSSatish Balay     isz[i] = isz_i;
375d5b485abSSatish Balay   }
3763a40ed3dSBarry Smith   PetscFunctionReturn(0);
377d5b485abSSatish Balay }
378d5b485abSSatish Balay /*
3792d4ee042Sprj-       MatIncreaseOverlap_MPIBAIJ_Receive - Process the received messages,
380d5b485abSSatish Balay          and return the output
381d5b485abSSatish Balay 
382d5b485abSSatish Balay          Input:
383d5b485abSSatish Balay            C    - the matrix
384d5b485abSSatish Balay            nrqr - no of messages being processed.
3852d4ee042Sprj-            rbuf - an array of pointers to the received requests
386d5b485abSSatish Balay 
387d5b485abSSatish Balay          Output:
388d5b485abSSatish Balay            xdata - array of messages to be sent back
389d5b485abSSatish Balay            isz1  - size of each message
390d5b485abSSatish Balay 
391a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
392d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
393a5b23f4aSJose E. Roman rather than all previous rows as it is now where a single large chunk of
394d5b485abSSatish Balay memory is used.
395d5b485abSSatish Balay 
396d5b485abSSatish Balay */
397b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
398d5b485abSSatish Balay {
399df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
400d5b485abSSatish Balay   Mat            A  = c->A,B = c->B;
401df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
402b24ad042SBarry Smith   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
403b24ad042SBarry Smith   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
404d2a212eaSBarry Smith   PetscInt       val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
405b24ad042SBarry Smith   PetscInt       *rbuf_i,kmax,rbuf_0;
406f1af5d2fSBarry Smith   PetscBT        xtable;
407d5b485abSSatish Balay 
4083a40ed3dSBarry Smith   PetscFunctionBegin;
409df36731bSBarry Smith   Mbs    = c->Mbs;
410899cda47SBarry Smith   rstart = c->rstartbs;
411899cda47SBarry Smith   cstart = c->cstartbs;
412d5b485abSSatish Balay   ai     = a->i;
413df36731bSBarry Smith   aj     = a->j;
414d5b485abSSatish Balay   bi     = b->i;
415df36731bSBarry Smith   bj     = b->j;
416d5b485abSSatish Balay   garray = c->garray;
417d5b485abSSatish Balay 
418d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
419d5b485abSSatish Balay     rbuf_i =  rbuf[i];
420d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
421d5b485abSSatish Balay     ct    += rbuf_0;
42226fbe8dcSKarl Rupp     for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j];
423d5b485abSSatish Balay   }
424d5b485abSSatish Balay 
425701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
426701b8814SBarry Smith   else        max1 = 1;
427d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
428175c2bc5SJunchao Zhang   if (nrqr) {
429*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(mem_estimate,&xdata[0]));
430d5b485abSSatish Balay     ++no_malloc;
431175c2bc5SJunchao Zhang   }
432*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTCreate(Mbs,&xtable));
433*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(isz1,nrqr));
434d5b485abSSatish Balay 
435d5b485abSSatish Balay   ct3 = 0;
436d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
437d5b485abSSatish Balay     rbuf_i =  rbuf[i];
438d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
439d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
440d5b485abSSatish Balay     ct2    =  ct1;
441d5b485abSSatish Balay     ct3   += ct1;
442d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
443*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTMemzero(Mbs,xtable));
444d5b485abSSatish Balay       oct2 = ct2;
445d5b485abSSatish Balay       kmax = rbuf_i[2*j];
446d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
447d5b485abSSatish Balay         row = rbuf_i[ct1];
448f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
449d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
450b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
451*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscMalloc1(new_estimate,&tmp));
452*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscArraycpy(tmp,xdata[0],mem_estimate));
453*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscFree(xdata[0]));
454d5b485abSSatish Balay             xdata[0]     = tmp;
455d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
45626fbe8dcSKarl Rupp             for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
457d5b485abSSatish Balay           }
458d5b485abSSatish Balay           xdata[i][ct2++] = row;
459d5b485abSSatish Balay           ct3++;
460d5b485abSSatish Balay         }
461d5b485abSSatish Balay       }
462d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
463d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
464d5b485abSSatish Balay         start = ai[row];
465d5b485abSSatish Balay         end   = ai[row+1];
466d5b485abSSatish Balay         for (l=start; l<end; l++) {
467df36731bSBarry Smith           val = aj[l] + cstart;
468f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
469d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
470b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
471*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscMalloc1(new_estimate,&tmp));
472*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscArraycpy(tmp,xdata[0],mem_estimate));
473*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscFree(xdata[0]));
474d5b485abSSatish Balay               xdata[0]     = tmp;
475d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
47626fbe8dcSKarl Rupp               for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
477d5b485abSSatish Balay             }
478d5b485abSSatish Balay             xdata[i][ct2++] = val;
479d5b485abSSatish Balay             ct3++;
480d5b485abSSatish Balay           }
481d5b485abSSatish Balay         }
482d5b485abSSatish Balay         start = bi[row];
483d5b485abSSatish Balay         end   = bi[row+1];
484d5b485abSSatish Balay         for (l=start; l<end; l++) {
485df36731bSBarry Smith           val = garray[bj[l]];
486f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
487d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
488b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
489*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscMalloc1(new_estimate,&tmp));
490*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscArraycpy(tmp,xdata[0],mem_estimate));
491*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscFree(xdata[0]));
492d5b485abSSatish Balay               xdata[0]     = tmp;
493d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
49426fbe8dcSKarl Rupp               for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
495d5b485abSSatish Balay             }
496d5b485abSSatish Balay             xdata[i][ct2++] = val;
497d5b485abSSatish Balay             ct3++;
498d5b485abSSatish Balay           }
499d5b485abSSatish Balay         }
500d5b485abSSatish Balay       }
501d5b485abSSatish Balay       /* Update the header*/
502d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
503d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
504d5b485abSSatish Balay     }
505d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
506175c2bc5SJunchao Zhang     if (i+1<nrqr) xdata[i+1] = xdata[i] + ct2;
507d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
508d5b485abSSatish Balay   }
509*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTDestroy(&xtable));
510*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(C,"Allocated %" PetscInt_FMT " bytes, required %" PetscInt_FMT ", no of mallocs = %" PetscInt_FMT "\n",mem_estimate,ct3,no_malloc));
5113a40ed3dSBarry Smith   PetscFunctionReturn(0);
512d5b485abSSatish Balay }
513d5b485abSSatish Balay 
5147dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
515d5b485abSSatish Balay {
51617df9f7cSHong Zhang   IS             *isrow_block,*iscol_block;
517cf4f527aSSatish Balay   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
518121971b2SHong Zhang   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,N=C->cmap->N,bs=C->rmap->bs;
5195dc98d6aSHong Zhang   Mat_SeqBAIJ    *subc;
5205c39f6d9SHong Zhang   Mat_SubSppt    *smat;
521a2feddc5SSatish Balay 
5223a40ed3dSBarry Smith   PetscFunctionBegin;
52329dcf524SDmitry Karpeev   /* The compression and expansion should be avoided. Doesn't point
52429dcf524SDmitry Karpeev      out errors, might change the indices, hence buggey */
525*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block));
526*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block));
527*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block));
528cf4f527aSSatish Balay 
529cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
5304698041cSHong Zhang   if (!C->cmap->N) nmax=20*1000000/sizeof(PetscInt);
5314698041cSHong Zhang   else nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt));
532329f5518SBarry Smith   if (!nmax) nmax = 1;
5335dc98d6aSHong Zhang 
5345dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
5354698041cSHong Zhang     nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */
536cf4f527aSSatish Balay 
537653e4784SBarry Smith     /* Make sure every processor loops through the nstages */
538*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C)));
5395dc98d6aSHong Zhang 
5404698041cSHong Zhang     /* Allocate memory to hold all the submatrices and dummy submatrices */
541*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(ismax+nstages,submat));
5424698041cSHong Zhang   } else { /* MAT_REUSE_MATRIX */
5434698041cSHong Zhang     if (ismax) {
5445dc98d6aSHong Zhang       subc = (Mat_SeqBAIJ*)((*submat)[0]->data);
5455dc98d6aSHong Zhang       smat   = subc->submatis1;
5464698041cSHong Zhang     } else { /* (*submat)[0] is a dummy matrix */
5475c39f6d9SHong Zhang       smat = (Mat_SubSppt*)(*submat)[0]->data;
5484698041cSHong Zhang     }
5492c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!smat,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat");
5505dc98d6aSHong Zhang     nstages = smat->nstages;
5515dc98d6aSHong Zhang   }
5525dc98d6aSHong Zhang 
553cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
554cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
5554e195584SJunchao Zhang     else if (pos >= ismax) max_no = 0;
556cf4f527aSSatish Balay     else                   max_no = ismax-pos;
557a42ab0d6SHong Zhang 
558*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,*submat+pos));
5594e195584SJunchao Zhang     if (!max_no) {
5604e195584SJunchao Zhang       if (scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */
5615c39f6d9SHong Zhang         smat = (Mat_SubSppt*)(*submat)[pos]->data;
5624698041cSHong Zhang         smat->nstages = nstages;
5634698041cSHong Zhang       }
5644e195584SJunchao Zhang       pos++; /* advance to next dummy matrix if any */
5654e195584SJunchao Zhang     } else pos += max_no;
566cf4f527aSSatish Balay   }
56736f4e84dSSatish Balay 
5685dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX && ismax) {
5695dc98d6aSHong Zhang     /* save nstages for reuse */
5705dc98d6aSHong Zhang     subc = (Mat_SeqBAIJ*)((*submat)[0]->data);
5715dc98d6aSHong Zhang     smat = subc->submatis1;
5725dc98d6aSHong Zhang     smat->nstages = nstages;
5735dc98d6aSHong Zhang   }
5745dc98d6aSHong Zhang 
57536f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
576*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&isrow_block[i]));
577*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&iscol_block[i]));
57836f4e84dSSatish Balay   }
579*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(isrow_block,iscol_block));
5803a40ed3dSBarry Smith   PetscFunctionReturn(0);
581a2feddc5SSatish Balay }
582a2feddc5SSatish Balay 
583233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
584e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
585233c2ff6SSatish Balay {
586e005ede5SBarry Smith   PetscInt       nGlobalNd = proc_gnode[size];
5874dc2109aSBarry Smith   PetscMPIInt    fproc;
588233c2ff6SSatish Balay 
589233c2ff6SSatish Balay   PetscFunctionBegin;
590*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc));
59123ce1328SBarry Smith   if (fproc > size) fproc = size;
592e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
593e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
594233c2ff6SSatish Balay     else                         fproc++;
595233c2ff6SSatish Balay   }
596e005ede5SBarry Smith   *rank = fproc;
597233c2ff6SSatish Balay   PetscFunctionReturn(0);
598233c2ff6SSatish Balay }
599233c2ff6SSatish Balay #endif
600233c2ff6SSatish Balay 
601a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
602b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
6037dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
60417df9f7cSHong Zhang {
60517df9f7cSHong Zhang   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
60617df9f7cSHong Zhang   Mat            A  = c->A;
60717df9f7cSHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc;
60817df9f7cSHong Zhang   const PetscInt **icol,**irow;
60917df9f7cSHong Zhang   PetscInt       *nrow,*ncol,start;
61017df9f7cSHong Zhang   PetscMPIInt    rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr;
611ddb3d6bcSHong Zhang   PetscInt       **sbuf1,**sbuf2,*sbuf2_i,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1;
61217df9f7cSHong Zhang   PetscInt       nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol;
61317df9f7cSHong Zhang   PetscInt       **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2;
61417df9f7cSHong Zhang   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
61517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
61617df9f7cSHong Zhang   PetscTable     *cmap,cmap_i=NULL,*rmap,rmap_i;
61717df9f7cSHong Zhang #else
61817df9f7cSHong Zhang   PetscInt       **cmap,*cmap_i=NULL,**rmap,*rmap_i;
61917df9f7cSHong Zhang #endif
62096cff833SHong Zhang   const PetscInt *irow_i,*icol_i;
62117df9f7cSHong Zhang   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
62217df9f7cSHong Zhang   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
62317df9f7cSHong Zhang   MPI_Request    *r_waits4,*s_waits3,*s_waits4;
62417df9f7cSHong Zhang   MPI_Comm       comm;
62564f5bb2dSSatish Balay   PetscScalar    **rbuf4,*rbuf4_i=NULL,**sbuf_aa,*vals,*mat_a=NULL,*imat_a=NULL,*sbuf_aa_i;
62617df9f7cSHong Zhang   PetscMPIInt    *onodes1,*olengths1,end;
62780d31651SHong Zhang   PetscInt       **row2proc,*row2proc_i,*imat_ilen,*imat_j,*imat_i;
6285c39f6d9SHong Zhang   Mat_SubSppt    *smat_i;
62917df9f7cSHong Zhang   PetscBool      *issorted,colflag,iscsorted=PETSC_TRUE;
63017df9f7cSHong Zhang   PetscInt       *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen;
631a42ab0d6SHong Zhang   PetscInt       bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs;
632a42ab0d6SHong Zhang   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
633a42ab0d6SHong Zhang   PetscInt       nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB;
634afb3cbd8SHong Zhang   PetscScalar    *vworkA=NULL,*vworkB=NULL,*a_a = a->a,*b_a = b->a;
635a42ab0d6SHong Zhang   PetscInt       cstart = c->cstartbs,*bmap = c->garray;
63680d31651SHong Zhang   PetscBool      *allrows,*allcolumns;
637a42ab0d6SHong Zhang 
63817df9f7cSHong Zhang   PetscFunctionBegin;
639*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)C,&comm));
64017df9f7cSHong Zhang   size = c->size;
64117df9f7cSHong Zhang   rank = c->rank;
64217df9f7cSHong Zhang 
643*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc5(ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax+1,&allcolumns,ismax,&allrows));
644*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc5(ismax,(PetscInt***)&irow,ismax,(PetscInt***)&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted));
64517df9f7cSHong Zhang 
64617df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
647*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSorted(iscol[i],&issorted[i]));
64880d31651SHong Zhang     if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */
649*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSorted(isrow[i],&issorted[i]));
65017df9f7cSHong Zhang 
651ec3c739cSHong Zhang     /* Check for special case: allcolumns */
652*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISIdentity(iscol[i],&colflag));
653*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(iscol[i],&ncol[i]));
6545dd0c08cSHong Zhang 
6551abc2f7bSHong Zhang     if (colflag && ncol[i] == c->Nbs) {
6561abc2f7bSHong Zhang       allcolumns[i] = PETSC_TRUE;
6571abc2f7bSHong Zhang       icol[i]       = NULL;
6581abc2f7bSHong Zhang     } else {
65917df9f7cSHong Zhang       allcolumns[i] = PETSC_FALSE;
660*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGetIndices(iscol[i],&icol[i]));
6611abc2f7bSHong Zhang     }
662ec3c739cSHong Zhang 
663ec3c739cSHong Zhang     /* Check for special case: allrows */
664*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISIdentity(isrow[i],&colflag));
665*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(isrow[i],&nrow[i]));
666ec3c739cSHong Zhang     if (colflag && nrow[i] == c->Mbs) {
667ec3c739cSHong Zhang       allrows[i] = PETSC_TRUE;
668ec3c739cSHong Zhang       irow[i]    = NULL;
669ec3c739cSHong Zhang     } else {
670ec3c739cSHong Zhang       allrows[i] = PETSC_FALSE;
671*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGetIndices(isrow[i],&irow[i]));
672ec3c739cSHong Zhang     }
67317df9f7cSHong Zhang   }
67417df9f7cSHong Zhang 
67517df9f7cSHong Zhang   if (scall == MAT_REUSE_MATRIX) {
67617df9f7cSHong Zhang     /* Assumes new rows are same length as the old rows */
67717df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
67817df9f7cSHong Zhang       subc = (Mat_SeqBAIJ*)(submats[i]->data);
6792c71b3e2SJacob Faibussowitsch       PetscCheckFalse(subc->mbs != nrow[i] || subc->nbs != ncol[i],PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
68017df9f7cSHong Zhang 
68117df9f7cSHong Zhang       /* Initial matrix as if empty */
682*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(subc->ilen,subc->mbs));
68317df9f7cSHong Zhang 
68417df9f7cSHong Zhang       /* Initial matrix as if empty */
68517df9f7cSHong Zhang       submats[i]->factortype = C->factortype;
68617df9f7cSHong Zhang 
68717df9f7cSHong Zhang       smat_i   = subc->submatis1;
68817df9f7cSHong Zhang 
68917df9f7cSHong Zhang       nrqs        = smat_i->nrqs;
69017df9f7cSHong Zhang       nrqr        = smat_i->nrqr;
69117df9f7cSHong Zhang       rbuf1       = smat_i->rbuf1;
69217df9f7cSHong Zhang       rbuf2       = smat_i->rbuf2;
69317df9f7cSHong Zhang       rbuf3       = smat_i->rbuf3;
69417df9f7cSHong Zhang       req_source2 = smat_i->req_source2;
69517df9f7cSHong Zhang 
69617df9f7cSHong Zhang       sbuf1     = smat_i->sbuf1;
69717df9f7cSHong Zhang       sbuf2     = smat_i->sbuf2;
69817df9f7cSHong Zhang       ptr       = smat_i->ptr;
69917df9f7cSHong Zhang       tmp       = smat_i->tmp;
70017df9f7cSHong Zhang       ctr       = smat_i->ctr;
70117df9f7cSHong Zhang 
70217df9f7cSHong Zhang       pa          = smat_i->pa;
70317df9f7cSHong Zhang       req_size    = smat_i->req_size;
70417df9f7cSHong Zhang       req_source1 = smat_i->req_source1;
70517df9f7cSHong Zhang 
70617df9f7cSHong Zhang       allcolumns[i] = smat_i->allcolumns;
707ec3c739cSHong Zhang       allrows[i]    = smat_i->allrows;
70817df9f7cSHong Zhang       row2proc[i]   = smat_i->row2proc;
70917df9f7cSHong Zhang       rmap[i]       = smat_i->rmap;
71017df9f7cSHong Zhang       cmap[i]       = smat_i->cmap;
71117df9f7cSHong Zhang     }
7124698041cSHong Zhang 
7134698041cSHong Zhang     if (!ismax){ /* Get dummy submatrices and retrieve struct submatis1 */
7142c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!submats[0],PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"submats are null, cannot reuse");
7155c39f6d9SHong Zhang       smat_i = (Mat_SubSppt*)submats[0]->data;
7164698041cSHong Zhang 
7174698041cSHong Zhang       nrqs        = smat_i->nrqs;
7184698041cSHong Zhang       nrqr        = smat_i->nrqr;
7194698041cSHong Zhang       rbuf1       = smat_i->rbuf1;
7204698041cSHong Zhang       rbuf2       = smat_i->rbuf2;
7214698041cSHong Zhang       rbuf3       = smat_i->rbuf3;
7224698041cSHong Zhang       req_source2 = smat_i->req_source2;
7234698041cSHong Zhang 
7244698041cSHong Zhang       sbuf1       = smat_i->sbuf1;
7254698041cSHong Zhang       sbuf2       = smat_i->sbuf2;
7264698041cSHong Zhang       ptr         = smat_i->ptr;
7274698041cSHong Zhang       tmp         = smat_i->tmp;
7284698041cSHong Zhang       ctr         = smat_i->ctr;
7294698041cSHong Zhang 
7304698041cSHong Zhang       pa          = smat_i->pa;
7314698041cSHong Zhang       req_size    = smat_i->req_size;
7324698041cSHong Zhang       req_source1 = smat_i->req_source1;
7334698041cSHong Zhang 
734c67fe082SHong Zhang       allcolumns[0] = PETSC_FALSE;
7354698041cSHong Zhang     }
73617df9f7cSHong Zhang   } else { /* scall == MAT_INITIAL_MATRIX */
73717df9f7cSHong Zhang     /* Get some new tags to keep the communication clean */
738*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetNewTag((PetscObject)C,&tag2));
739*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetNewTag((PetscObject)C,&tag3));
74017df9f7cSHong Zhang 
74117df9f7cSHong Zhang     /* evaluate communication - mesg to who, length of mesg, and buffer space
74217df9f7cSHong Zhang      required. Based on this, buffers are allocated, and data copied into them*/
743*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4));   /* mesg size, initialize work vectors */
74417df9f7cSHong Zhang 
74517df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
74617df9f7cSHong Zhang       jmax   = nrow[i];
74717df9f7cSHong Zhang       irow_i = irow[i];
74817df9f7cSHong Zhang 
749*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(jmax,&row2proc_i));
75017df9f7cSHong Zhang       row2proc[i] = row2proc_i;
75117df9f7cSHong Zhang 
75217df9f7cSHong Zhang       if (issorted[i]) proc = 0;
75317df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
75417df9f7cSHong Zhang         if (!issorted[i]) proc = 0;
755ec3c739cSHong Zhang         if (allrows[i]) row = j;
756ec3c739cSHong Zhang         else row = irow_i[j];
757ec3c739cSHong Zhang 
758a42ab0d6SHong Zhang         while (row >= c->rangebs[proc+1]) proc++;
75917df9f7cSHong Zhang         w4[proc]++;
76017df9f7cSHong Zhang         row2proc_i[j] = proc; /* map row index to proc */
76117df9f7cSHong Zhang       }
76217df9f7cSHong Zhang       for (j=0; j<size; j++) {
76317df9f7cSHong Zhang         if (w4[j]) { w1[j] += w4[j];  w3[j]++; w4[j] = 0;}
76417df9f7cSHong Zhang       }
76517df9f7cSHong Zhang     }
76617df9f7cSHong Zhang 
76717df9f7cSHong Zhang     nrqs     = 0;              /* no of outgoing messages */
76817df9f7cSHong Zhang     msz      = 0;              /* total mesg length (for all procs) */
76917df9f7cSHong Zhang     w1[rank] = 0;              /* no mesg sent to self */
77017df9f7cSHong Zhang     w3[rank] = 0;
77117df9f7cSHong Zhang     for (i=0; i<size; i++) {
77217df9f7cSHong Zhang       if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
77317df9f7cSHong Zhang     }
774*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqs,&pa)); /*(proc -array)*/
77517df9f7cSHong Zhang     for (i=0,j=0; i<size; i++) {
77617df9f7cSHong Zhang       if (w1[i]) { pa[j] = i; j++; }
77717df9f7cSHong Zhang     }
77817df9f7cSHong Zhang 
77917df9f7cSHong Zhang     /* Each message would have a header = 1 + 2*(no of IS) + data */
78017df9f7cSHong Zhang     for (i=0; i<nrqs; i++) {
78117df9f7cSHong Zhang       j      = pa[i];
78217df9f7cSHong Zhang       w1[j] += w2[j] + 2* w3[j];
78317df9f7cSHong Zhang       msz   += w1[j];
78417df9f7cSHong Zhang     }
785*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(0,"Number of outgoing messages %" PetscInt_FMT " Total message length %" PetscInt_FMT "\n",nrqs,msz));
78617df9f7cSHong Zhang 
78717df9f7cSHong Zhang     /* Determine the number of messages to expect, their lengths, from from-ids */
788*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscGatherNumberOfMessages(comm,w2,w1,&nrqr));
789*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1));
79017df9f7cSHong Zhang 
79117df9f7cSHong Zhang     /* Now post the Irecvs corresponding to these messages */
792*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetNewTag((PetscObject)C,&tag0));
793*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1));
79417df9f7cSHong Zhang 
79517df9f7cSHong Zhang     /* Allocate Memory for outgoing messages */
796*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr));
797*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(sbuf1,size));
798*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(ptr,size));
79917df9f7cSHong Zhang 
80017df9f7cSHong Zhang     {
80117df9f7cSHong Zhang       PetscInt *iptr = tmp;
80217df9f7cSHong Zhang       k    = 0;
80317df9f7cSHong Zhang       for (i=0; i<nrqs; i++) {
80417df9f7cSHong Zhang         j        = pa[i];
80517df9f7cSHong Zhang         iptr    += k;
80617df9f7cSHong Zhang         sbuf1[j] = iptr;
80717df9f7cSHong Zhang         k        = w1[j];
80817df9f7cSHong Zhang       }
80917df9f7cSHong Zhang     }
81017df9f7cSHong Zhang 
81117df9f7cSHong Zhang     /* Form the outgoing messages. Initialize the header space */
81217df9f7cSHong Zhang     for (i=0; i<nrqs; i++) {
81317df9f7cSHong Zhang       j           = pa[i];
81417df9f7cSHong Zhang       sbuf1[j][0] = 0;
815*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(sbuf1[j]+1,2*w3[j]));
81617df9f7cSHong Zhang       ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
81717df9f7cSHong Zhang     }
81817df9f7cSHong Zhang 
81917df9f7cSHong Zhang     /* Parse the isrow and copy data into outbuf */
82017df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
82117df9f7cSHong Zhang       row2proc_i = row2proc[i];
822*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctr,size));
82317df9f7cSHong Zhang       irow_i = irow[i];
82417df9f7cSHong Zhang       jmax   = nrow[i];
82517df9f7cSHong Zhang       for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
82617df9f7cSHong Zhang         proc = row2proc_i[j];
827ec3c739cSHong Zhang         if (allrows[i]) row = j;
828ec3c739cSHong Zhang         else row = irow_i[j];
829ec3c739cSHong Zhang 
83017df9f7cSHong Zhang         if (proc != rank) { /* copy to the outgoing buf*/
83117df9f7cSHong Zhang           ctr[proc]++;
832ec3c739cSHong Zhang           *ptr[proc] = row;
83317df9f7cSHong Zhang           ptr[proc]++;
83417df9f7cSHong Zhang         }
83517df9f7cSHong Zhang       }
83617df9f7cSHong Zhang       /* Update the headers for the current IS */
83717df9f7cSHong Zhang       for (j=0; j<size; j++) { /* Can Optimise this loop too */
83817df9f7cSHong Zhang         if ((ctr_j = ctr[j])) {
83917df9f7cSHong Zhang           sbuf1_j        = sbuf1[j];
84017df9f7cSHong Zhang           k              = ++sbuf1_j[0];
84117df9f7cSHong Zhang           sbuf1_j[2*k]   = ctr_j;
84217df9f7cSHong Zhang           sbuf1_j[2*k-1] = i;
84317df9f7cSHong Zhang         }
84417df9f7cSHong Zhang       }
84517df9f7cSHong Zhang     }
84617df9f7cSHong Zhang 
84717df9f7cSHong Zhang     /*  Now  post the sends */
848*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqs,&s_waits1));
84917df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
85017df9f7cSHong Zhang       j    = pa[i];
851*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i));
85217df9f7cSHong Zhang     }
85317df9f7cSHong Zhang 
85417df9f7cSHong Zhang     /* Post Receives to capture the buffer size */
855*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqs,&r_waits2));
856*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc3(nrqs,&req_source2,nrqs,&rbuf2,nrqs,&rbuf3));
857175c2bc5SJunchao Zhang     if (nrqs) rbuf2[0] = tmp + msz;
85817df9f7cSHong Zhang     for (i=1; i<nrqs; ++i) {
85917df9f7cSHong Zhang       rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
86017df9f7cSHong Zhang     }
86117df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
86217df9f7cSHong Zhang       j    = pa[i];
863*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i));
86417df9f7cSHong Zhang     }
86517df9f7cSHong Zhang 
86617df9f7cSHong Zhang     /* Send to other procs the buf size they should allocate */
86717df9f7cSHong Zhang     /* Receive messages*/
868*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqr,&s_waits2));
869*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1));
87017df9f7cSHong Zhang 
871*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrqr,r_waits1,MPI_STATUSES_IGNORE));
87217df9f7cSHong Zhang     for (i=0; i<nrqr; ++i) {
87317df9f7cSHong Zhang       req_size[i] = 0;
87417df9f7cSHong Zhang       rbuf1_i        = rbuf1[i];
87517df9f7cSHong Zhang       start          = 2*rbuf1_i[0] + 1;
876175c2bc5SJunchao Zhang       end            = olengths1[i];
877*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(end,&sbuf2[i]));
87817df9f7cSHong Zhang       sbuf2_i        = sbuf2[i];
87917df9f7cSHong Zhang       for (j=start; j<end; j++) {
880ddb3d6bcSHong Zhang         row             = rbuf1_i[j] - rstart;
881ddb3d6bcSHong Zhang         ncols           = a_i[row+1] - a_i[row] + b_i[row+1] - b_i[row];
88217df9f7cSHong Zhang         sbuf2_i[j]      = ncols;
88317df9f7cSHong Zhang         req_size[i] += ncols;
88417df9f7cSHong Zhang       }
885175c2bc5SJunchao Zhang       req_source1[i] = onodes1[i];
88617df9f7cSHong Zhang       /* form the header */
88717df9f7cSHong Zhang       sbuf2_i[0] = req_size[i];
88817df9f7cSHong Zhang       for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
88917df9f7cSHong Zhang 
890*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i));
89117df9f7cSHong Zhang     }
892ddb3d6bcSHong Zhang 
893*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(onodes1));
894*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(olengths1));
895175c2bc5SJunchao Zhang 
896*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(r_waits1));
897*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree4(w1,w2,w3,w4));
89817df9f7cSHong Zhang 
89917df9f7cSHong Zhang     /* Receive messages*/
900*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqs,&r_waits3));
90117df9f7cSHong Zhang 
902*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrqs,r_waits2,MPI_STATUSES_IGNORE));
90317df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
904*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(rbuf2[i][0],&rbuf3[i]));
905175c2bc5SJunchao Zhang       req_source2[i] = pa[i];
906*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i));
90717df9f7cSHong Zhang     }
908*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(r_waits2));
90917df9f7cSHong Zhang 
91017df9f7cSHong Zhang     /* Wait on sends1 and sends2 */
911*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrqs,s_waits1,MPI_STATUSES_IGNORE));
912*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrqr,s_waits2,MPI_STATUSES_IGNORE));
913*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(s_waits1));
914*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(s_waits2));
91517df9f7cSHong Zhang 
91617df9f7cSHong Zhang     /* Now allocate sending buffers for a->j, and send them off */
917*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqr,&sbuf_aj));
91817df9f7cSHong Zhang     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
919*5f80ce2aSJacob Faibussowitsch     if (nrqr) CHKERRQ(PetscMalloc1(j,&sbuf_aj[0]));
92017df9f7cSHong Zhang     for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
92117df9f7cSHong Zhang 
922*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqr,&s_waits3));
92317df9f7cSHong Zhang     {
92417df9f7cSHong Zhang 
92517df9f7cSHong Zhang       for (i=0; i<nrqr; i++) {
92617df9f7cSHong Zhang         rbuf1_i   = rbuf1[i];
92717df9f7cSHong Zhang         sbuf_aj_i = sbuf_aj[i];
92817df9f7cSHong Zhang         ct1       = 2*rbuf1_i[0] + 1;
92917df9f7cSHong Zhang         ct2       = 0;
93017df9f7cSHong Zhang         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
93117df9f7cSHong Zhang           kmax = rbuf1[i][2*j];
93217df9f7cSHong Zhang           for (k=0; k<kmax; k++,ct1++) {
93317df9f7cSHong Zhang             row    = rbuf1_i[ct1] - rstart;
93417df9f7cSHong Zhang             nzA    = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row];
93517df9f7cSHong Zhang             ncols  = nzA + nzB;
93617df9f7cSHong Zhang             cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
93717df9f7cSHong Zhang 
93817df9f7cSHong Zhang             /* load the column indices for this row into cols */
93917df9f7cSHong Zhang             cols = sbuf_aj_i + ct2;
94017df9f7cSHong Zhang             for (l=0; l<nzB; l++) {
941a42ab0d6SHong Zhang               if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
942a42ab0d6SHong Zhang               else break;
94317df9f7cSHong Zhang             }
944a42ab0d6SHong Zhang             imark = l;
945a42ab0d6SHong Zhang             for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];}
946a42ab0d6SHong Zhang             for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
94717df9f7cSHong Zhang             ct2 += ncols;
94817df9f7cSHong Zhang           }
94917df9f7cSHong Zhang         }
950*5f80ce2aSJacob Faibussowitsch         CHKERRMPI(MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i));
95117df9f7cSHong Zhang       }
95217df9f7cSHong Zhang     }
95317df9f7cSHong Zhang 
95417df9f7cSHong Zhang     /* create col map: global col of C -> local col of submatrices */
95517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
95617df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
95717df9f7cSHong Zhang       if (!allcolumns[i]) {
958*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscTableCreate(ncol[i],c->Nbs,&cmap[i]));
95917df9f7cSHong Zhang 
96017df9f7cSHong Zhang         jmax   = ncol[i];
96117df9f7cSHong Zhang         icol_i = icol[i];
96217df9f7cSHong Zhang         cmap_i = cmap[i];
96317df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
964*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES));
96517df9f7cSHong Zhang         }
96617df9f7cSHong Zhang       } else cmap[i] = NULL;
96717df9f7cSHong Zhang     }
96817df9f7cSHong Zhang #else
96917df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
97017df9f7cSHong Zhang       if (!allcolumns[i]) {
971*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscCalloc1(c->Nbs,&cmap[i]));
97217df9f7cSHong Zhang         jmax   = ncol[i];
97317df9f7cSHong Zhang         icol_i = icol[i];
97417df9f7cSHong Zhang         cmap_i = cmap[i];
975a42ab0d6SHong Zhang         for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1;
97617df9f7cSHong Zhang       } else cmap[i] = NULL;
97717df9f7cSHong Zhang     }
97817df9f7cSHong Zhang #endif
97917df9f7cSHong Zhang 
98017df9f7cSHong Zhang     /* Create lens which is required for MatCreate... */
98117df9f7cSHong Zhang     for (i=0,j=0; i<ismax; i++) j += nrow[i];
982*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(ismax,&lens));
98317df9f7cSHong Zhang 
98417df9f7cSHong Zhang     if (ismax) {
985*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc1(j,&lens[0]));
98617df9f7cSHong Zhang     }
98717df9f7cSHong Zhang     for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];
98817df9f7cSHong Zhang 
98917df9f7cSHong Zhang     /* Update lens from local data */
99017df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
99117df9f7cSHong Zhang       row2proc_i = row2proc[i];
99217df9f7cSHong Zhang       jmax = nrow[i];
99317df9f7cSHong Zhang       if (!allcolumns[i]) cmap_i = cmap[i];
99417df9f7cSHong Zhang       irow_i = irow[i];
99517df9f7cSHong Zhang       lens_i = lens[i];
99617df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
997ec3c739cSHong Zhang         if (allrows[i]) row = j;
998ec3c739cSHong Zhang         else row = irow_i[j]; /* global blocked row of C */
999ec3c739cSHong Zhang 
100017df9f7cSHong Zhang         proc = row2proc_i[j];
100117df9f7cSHong Zhang         if (proc == rank) {
1002a42ab0d6SHong Zhang           /* Get indices from matA and then from matB */
100317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1004a42ab0d6SHong Zhang           PetscInt   tt;
1005a42ab0d6SHong Zhang #endif
1006a42ab0d6SHong Zhang           row    = row - rstart;
1007a42ab0d6SHong Zhang           nzA    = a_i[row+1] - a_i[row];
1008a42ab0d6SHong Zhang           nzB    = b_i[row+1] - b_i[row];
1009a42ab0d6SHong Zhang           cworkA =  a_j + a_i[row];
1010a42ab0d6SHong Zhang           cworkB = b_j + b_i[row];
1011a42ab0d6SHong Zhang 
1012a42ab0d6SHong Zhang           if (!allcolumns[i]) {
1013a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1014a42ab0d6SHong Zhang             for (k=0; k<nzA; k++) {
1015*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt));
1016a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1017a42ab0d6SHong Zhang             }
1018a42ab0d6SHong Zhang             for (k=0; k<nzB; k++) {
1019*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt));
1020a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1021a42ab0d6SHong Zhang             }
1022a42ab0d6SHong Zhang 
1023a42ab0d6SHong Zhang #else
1024a42ab0d6SHong Zhang             for (k=0; k<nzA; k++) {
1025a42ab0d6SHong Zhang               if (cmap_i[cstart + cworkA[k]]) lens_i[j]++;
1026a42ab0d6SHong Zhang             }
1027a42ab0d6SHong Zhang             for (k=0; k<nzB; k++) {
1028a42ab0d6SHong Zhang               if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++;
1029a42ab0d6SHong Zhang             }
1030a42ab0d6SHong Zhang #endif
1031a42ab0d6SHong Zhang           } else { /* allcolumns */
1032a42ab0d6SHong Zhang             lens_i[j] = nzA + nzB;
1033a42ab0d6SHong Zhang           }
103417df9f7cSHong Zhang         }
103517df9f7cSHong Zhang       }
103617df9f7cSHong Zhang     }
103717df9f7cSHong Zhang 
103817df9f7cSHong Zhang     /* Create row map: global row of C -> local row of submatrices */
103917df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
1040059f6b74SHong Zhang       if (!allrows[i]) {
1041059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE)
1042*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscTableCreate(nrow[i],c->Mbs,&rmap[i]));
104317df9f7cSHong Zhang         irow_i = irow[i];
104417df9f7cSHong Zhang         jmax   = nrow[i];
104517df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
1046ec3c739cSHong Zhang           if (allrows[i]) {
1047*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscTableAdd(rmap[i],j+1,j+1,INSERT_VALUES));
1048ec3c739cSHong Zhang           } else {
1049*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES));
105017df9f7cSHong Zhang           }
105117df9f7cSHong Zhang         }
105217df9f7cSHong Zhang #else
1053*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscCalloc1(c->Mbs,&rmap[i]));
105417df9f7cSHong Zhang         rmap_i = rmap[i];
105517df9f7cSHong Zhang         irow_i = irow[i];
105617df9f7cSHong Zhang         jmax   = nrow[i];
105717df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
1058ec3c739cSHong Zhang           if (allrows[i]) rmap_i[j] = j;
1059ec3c739cSHong Zhang           else rmap_i[irow_i[j]] = j;
106017df9f7cSHong Zhang         }
106117df9f7cSHong Zhang #endif
1062059f6b74SHong Zhang       } else rmap[i] = NULL;
1063059f6b74SHong Zhang     }
106417df9f7cSHong Zhang 
106517df9f7cSHong Zhang     /* Update lens from offproc data */
106617df9f7cSHong Zhang     {
106717df9f7cSHong Zhang       PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i;
106817df9f7cSHong Zhang 
1069*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Waitall(nrqs,r_waits3,MPI_STATUSES_IGNORE));
107017df9f7cSHong Zhang       for (tmp2=0; tmp2<nrqs; tmp2++) {
107117df9f7cSHong Zhang         sbuf1_i = sbuf1[pa[tmp2]];
107217df9f7cSHong Zhang         jmax    = sbuf1_i[0];
107317df9f7cSHong Zhang         ct1     = 2*jmax+1;
107417df9f7cSHong Zhang         ct2     = 0;
107517df9f7cSHong Zhang         rbuf2_i = rbuf2[tmp2];
107617df9f7cSHong Zhang         rbuf3_i = rbuf3[tmp2];
107717df9f7cSHong Zhang         for (j=1; j<=jmax; j++) {
107817df9f7cSHong Zhang           is_no  = sbuf1_i[2*j-1];
107917df9f7cSHong Zhang           max1   = sbuf1_i[2*j];
108017df9f7cSHong Zhang           lens_i = lens[is_no];
108117df9f7cSHong Zhang           if (!allcolumns[is_no]) cmap_i = cmap[is_no];
108217df9f7cSHong Zhang           rmap_i = rmap[is_no];
108317df9f7cSHong Zhang           for (k=0; k<max1; k++,ct1++) {
1084059f6b74SHong Zhang             if (allrows[is_no]) {
1085059f6b74SHong Zhang               row = sbuf1_i[ct1];
1086059f6b74SHong Zhang             } else {
108717df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1088*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row));
108917df9f7cSHong Zhang               row--;
10902c71b3e2SJacob Faibussowitsch               PetscCheckFalse(row < 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
109117df9f7cSHong Zhang #else
109217df9f7cSHong Zhang               row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
109317df9f7cSHong Zhang #endif
1094059f6b74SHong Zhang             }
109517df9f7cSHong Zhang             max2 = rbuf2_i[ct1];
109617df9f7cSHong Zhang             for (l=0; l<max2; l++,ct2++) {
109717df9f7cSHong Zhang               if (!allcolumns[is_no]) {
109817df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1099*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol));
110017df9f7cSHong Zhang #else
110117df9f7cSHong Zhang                 tcol = cmap_i[rbuf3_i[ct2]];
110217df9f7cSHong Zhang #endif
110317df9f7cSHong Zhang                 if (tcol) lens_i[row]++;
110417df9f7cSHong Zhang               } else { /* allcolumns */
110517df9f7cSHong Zhang                 lens_i[row]++; /* lens_i[row] += max2 ? */
110617df9f7cSHong Zhang               }
110717df9f7cSHong Zhang             }
110817df9f7cSHong Zhang           }
110917df9f7cSHong Zhang         }
111017df9f7cSHong Zhang       }
111117df9f7cSHong Zhang     }
1112*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(r_waits3));
1113*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrqr,s_waits3,MPI_STATUSES_IGNORE));
1114*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(s_waits3));
111517df9f7cSHong Zhang 
111617df9f7cSHong Zhang     /* Create the submatrices */
111717df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
1118a42ab0d6SHong Zhang       PetscInt bs_tmp;
1119a42ab0d6SHong Zhang       if (ijonly) bs_tmp = 1;
1120a42ab0d6SHong Zhang       else        bs_tmp = bs;
112117df9f7cSHong Zhang 
1122*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreate(PETSC_COMM_SELF,submats+i));
1123*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE));
112417df9f7cSHong Zhang 
1125*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetType(submats[i],((PetscObject)A)->type_name));
1126*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]));
1127*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i])); /* this subroutine is used by SBAIJ routines */
112817df9f7cSHong Zhang 
11295c39f6d9SHong Zhang       /* create struct Mat_SubSppt and attached it to submat */
1130*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscNew(&smat_i));
113117df9f7cSHong Zhang       subc = (Mat_SeqBAIJ*)submats[i]->data;
113217df9f7cSHong Zhang       subc->submatis1 = smat_i;
113317df9f7cSHong Zhang 
113417df9f7cSHong Zhang       smat_i->destroy          = submats[i]->ops->destroy;
1135f68bb481SHong Zhang       submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ;
113617df9f7cSHong Zhang       submats[i]->factortype   = C->factortype;
113717df9f7cSHong Zhang 
113817df9f7cSHong Zhang       smat_i->id          = i;
113917df9f7cSHong Zhang       smat_i->nrqs        = nrqs;
114017df9f7cSHong Zhang       smat_i->nrqr        = nrqr;
114117df9f7cSHong Zhang       smat_i->rbuf1       = rbuf1;
114217df9f7cSHong Zhang       smat_i->rbuf2       = rbuf2;
114317df9f7cSHong Zhang       smat_i->rbuf3       = rbuf3;
114417df9f7cSHong Zhang       smat_i->sbuf2       = sbuf2;
114517df9f7cSHong Zhang       smat_i->req_source2 = req_source2;
114617df9f7cSHong Zhang 
114717df9f7cSHong Zhang       smat_i->sbuf1       = sbuf1;
114817df9f7cSHong Zhang       smat_i->ptr         = ptr;
114917df9f7cSHong Zhang       smat_i->tmp         = tmp;
115017df9f7cSHong Zhang       smat_i->ctr         = ctr;
115117df9f7cSHong Zhang 
115217df9f7cSHong Zhang       smat_i->pa           = pa;
115317df9f7cSHong Zhang       smat_i->req_size     = req_size;
115417df9f7cSHong Zhang       smat_i->req_source1  = req_source1;
115517df9f7cSHong Zhang 
115617df9f7cSHong Zhang       smat_i->allcolumns  = allcolumns[i];
1157ec3c739cSHong Zhang       smat_i->allrows     = allrows[i];
115817df9f7cSHong Zhang       smat_i->singleis    = PETSC_FALSE;
115917df9f7cSHong Zhang       smat_i->row2proc    = row2proc[i];
116017df9f7cSHong Zhang       smat_i->rmap        = rmap[i];
116117df9f7cSHong Zhang       smat_i->cmap        = cmap[i];
116217df9f7cSHong Zhang     }
116317df9f7cSHong Zhang 
11644698041cSHong Zhang     if (!ismax) { /* Create dummy submats[0] for reuse struct subc */
1165*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreate(PETSC_COMM_SELF,&submats[0]));
1166*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(submats[0],0,0,PETSC_DETERMINE,PETSC_DETERMINE));
1167*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetType(submats[0],MATDUMMY));
11684698041cSHong Zhang 
11695c39f6d9SHong Zhang       /* create struct Mat_SubSppt and attached it to submat */
1170*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscNewLog(submats[0],&smat_i));
11714698041cSHong Zhang       submats[0]->data = (void*)smat_i;
11724698041cSHong Zhang 
11734698041cSHong Zhang       smat_i->destroy          = submats[0]->ops->destroy;
1174f68bb481SHong Zhang       submats[0]->ops->destroy = MatDestroySubMatrix_Dummy;
11754698041cSHong Zhang       submats[0]->factortype   = C->factortype;
11764698041cSHong Zhang 
1177006cef31SHong Zhang       smat_i->id          = 0;
11784698041cSHong Zhang       smat_i->nrqs        = nrqs;
11794698041cSHong Zhang       smat_i->nrqr        = nrqr;
11804698041cSHong Zhang       smat_i->rbuf1       = rbuf1;
11814698041cSHong Zhang       smat_i->rbuf2       = rbuf2;
11824698041cSHong Zhang       smat_i->rbuf3       = rbuf3;
11834698041cSHong Zhang       smat_i->sbuf2       = sbuf2;
11844698041cSHong Zhang       smat_i->req_source2 = req_source2;
11854698041cSHong Zhang 
11864698041cSHong Zhang       smat_i->sbuf1       = sbuf1;
11874698041cSHong Zhang       smat_i->ptr         = ptr;
11884698041cSHong Zhang       smat_i->tmp         = tmp;
11894698041cSHong Zhang       smat_i->ctr         = ctr;
11904698041cSHong Zhang 
11914698041cSHong Zhang       smat_i->pa           = pa;
11924698041cSHong Zhang       smat_i->req_size     = req_size;
11934698041cSHong Zhang       smat_i->req_source1  = req_source1;
11944698041cSHong Zhang 
1195c67fe082SHong Zhang       smat_i->allcolumns  = PETSC_FALSE;
11964698041cSHong Zhang       smat_i->singleis    = PETSC_FALSE;
11974698041cSHong Zhang       smat_i->row2proc    = NULL;
11984698041cSHong Zhang       smat_i->rmap        = NULL;
11994698041cSHong Zhang       smat_i->cmap        = NULL;
12004698041cSHong Zhang     }
12014698041cSHong Zhang 
1202*5f80ce2aSJacob Faibussowitsch     if (ismax) CHKERRQ(PetscFree(lens[0]));
1203*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(lens));
1204175c2bc5SJunchao Zhang     if (sbuf_aj) {
1205*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(sbuf_aj[0]));
1206*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(sbuf_aj));
1207175c2bc5SJunchao Zhang     }
120817df9f7cSHong Zhang 
120917df9f7cSHong Zhang   } /* endof scall == MAT_INITIAL_MATRIX */
121017df9f7cSHong Zhang 
121117df9f7cSHong Zhang   /* Post recv matrix values */
1212bbc89d27SHong Zhang   if (!ijonly) {
1213*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetNewTag((PetscObject)C,&tag4));
1214*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqs,&rbuf4));
1215*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqs,&r_waits4));
121617df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
1217*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]));
1218*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i));
121917df9f7cSHong Zhang     }
122017df9f7cSHong Zhang 
122117df9f7cSHong Zhang     /* Allocate sending buffers for a->a, and send them off */
1222*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqr,&sbuf_aa));
122317df9f7cSHong Zhang     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
12249e686edcSHong Zhang 
1225*5f80ce2aSJacob Faibussowitsch     if (nrqr) CHKERRQ(PetscMalloc1(j*bs2,&sbuf_aa[0]));
1226a42ab0d6SHong Zhang     for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
122717df9f7cSHong Zhang 
1228*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrqr,&s_waits4));
122917df9f7cSHong Zhang 
123017df9f7cSHong Zhang     for (i=0; i<nrqr; i++) {
123117df9f7cSHong Zhang       rbuf1_i   = rbuf1[i];
123217df9f7cSHong Zhang       sbuf_aa_i = sbuf_aa[i];
123317df9f7cSHong Zhang       ct1       = 2*rbuf1_i[0]+1;
123417df9f7cSHong Zhang       ct2       = 0;
123517df9f7cSHong Zhang       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
123617df9f7cSHong Zhang         kmax = rbuf1_i[2*j];
123717df9f7cSHong Zhang         for (k=0; k<kmax; k++,ct1++) {
123817df9f7cSHong Zhang           row    = rbuf1_i[ct1] - rstart;
1239a42ab0d6SHong Zhang           nzA    = a_i[row+1] - a_i[row];
1240a42ab0d6SHong Zhang           nzB    = b_i[row+1] - b_i[row];
124117df9f7cSHong Zhang           ncols  = nzA + nzB;
124217df9f7cSHong Zhang           cworkB = b_j + b_i[row];
1243a42ab0d6SHong Zhang           vworkA = a_a + a_i[row]*bs2;
1244a42ab0d6SHong Zhang           vworkB = b_a + b_i[row]*bs2;
124517df9f7cSHong Zhang 
124617df9f7cSHong Zhang           /* load the column values for this row into vals*/
1247a42ab0d6SHong Zhang           vals = sbuf_aa_i+ct2*bs2;
1248a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1249a42ab0d6SHong Zhang             if ((bmap[cworkB[l]]) < cstart) {
1250*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscArraycpy(vals+l*bs2,vworkB+l*bs2,bs2));
1251a42ab0d6SHong Zhang             } else break;
1252a42ab0d6SHong Zhang           }
1253a42ab0d6SHong Zhang           imark = l;
1254a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1255*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscArraycpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2));
1256a42ab0d6SHong Zhang           }
1257a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1258*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscArraycpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2));
1259a42ab0d6SHong Zhang           }
126017df9f7cSHong Zhang 
126117df9f7cSHong Zhang           ct2 += ncols;
126217df9f7cSHong Zhang         }
126317df9f7cSHong Zhang       }
1264*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i));
126517df9f7cSHong Zhang     }
1266bbc89d27SHong Zhang   }
126717df9f7cSHong Zhang 
126817df9f7cSHong Zhang   /* Assemble the matrices */
126917df9f7cSHong Zhang   /* First assemble the local rows */
127017df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
127117df9f7cSHong Zhang     row2proc_i = row2proc[i];
127217df9f7cSHong Zhang     subc      = (Mat_SeqBAIJ*)submats[i]->data;
127317df9f7cSHong Zhang     imat_ilen = subc->ilen;
127417df9f7cSHong Zhang     imat_j    = subc->j;
127517df9f7cSHong Zhang     imat_i    = subc->i;
127617df9f7cSHong Zhang     imat_a    = subc->a;
127717df9f7cSHong Zhang 
127817df9f7cSHong Zhang     if (!allcolumns[i]) cmap_i = cmap[i];
127917df9f7cSHong Zhang     rmap_i = rmap[i];
128017df9f7cSHong Zhang     irow_i = irow[i];
128117df9f7cSHong Zhang     jmax   = nrow[i];
128217df9f7cSHong Zhang     for (j=0; j<jmax; j++) {
1283ec3c739cSHong Zhang       if (allrows[i]) row = j;
1284ec3c739cSHong Zhang       else row  = irow_i[j];
128517df9f7cSHong Zhang       proc = row2proc_i[j];
1286a42ab0d6SHong Zhang 
128717df9f7cSHong Zhang       if (proc == rank) {
1288a42ab0d6SHong Zhang 
1289a42ab0d6SHong Zhang         row    = row - rstart;
1290a42ab0d6SHong Zhang         nzA    = a_i[row+1] - a_i[row];
1291a42ab0d6SHong Zhang         nzB    = b_i[row+1] - b_i[row];
1292a42ab0d6SHong Zhang         cworkA = a_j + a_i[row];
1293a42ab0d6SHong Zhang         cworkB = b_j + b_i[row];
1294a42ab0d6SHong Zhang         if (!ijonly) {
1295a42ab0d6SHong Zhang           vworkA = a_a + a_i[row]*bs2;
1296a42ab0d6SHong Zhang           vworkB = b_a + b_i[row]*bs2;
1297a42ab0d6SHong Zhang         }
1298059f6b74SHong Zhang 
1299059f6b74SHong Zhang         if (allrows[i]) {
1300059f6b74SHong Zhang           row = row+rstart;
1301059f6b74SHong Zhang         } else {
1302a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1303*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscTableFind(rmap_i,row+rstart+1,&row));
1304a42ab0d6SHong Zhang           row--;
1305121971b2SHong Zhang 
13062c71b3e2SJacob Faibussowitsch           PetscCheckFalse(row < 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1307a42ab0d6SHong Zhang #else
1308a42ab0d6SHong Zhang           row = rmap_i[row + rstart];
1309a42ab0d6SHong Zhang #endif
1310059f6b74SHong Zhang         }
1311a42ab0d6SHong Zhang         mat_i = imat_i[row];
1312a42ab0d6SHong Zhang         if (!ijonly) mat_a = imat_a + mat_i*bs2;
1313a42ab0d6SHong Zhang         mat_j    = imat_j + mat_i;
131480d31651SHong Zhang         ilen = imat_ilen[row];
1315a42ab0d6SHong Zhang 
1316a42ab0d6SHong Zhang         /* load the column indices for this row into cols*/
1317a42ab0d6SHong Zhang         if (!allcolumns[i]) {
1318a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1319a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1320a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1321*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscTableFind(cmap_i,ctmp+1,&tcol));
1322a42ab0d6SHong Zhang               if (tcol) {
1323a42ab0d6SHong Zhang #else
1324a42ab0d6SHong Zhang               if ((tcol = cmap_i[ctmp])) {
1325a42ab0d6SHong Zhang #endif
1326a42ab0d6SHong Zhang                 *mat_j++ = tcol - 1;
1327*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(PetscArraycpy(mat_a,vworkB+l*bs2,bs2));
1328a42ab0d6SHong Zhang                 mat_a   += bs2;
132980d31651SHong Zhang                 ilen++;
1330a42ab0d6SHong Zhang               }
1331a42ab0d6SHong Zhang             } else break;
1332a42ab0d6SHong Zhang           }
1333a42ab0d6SHong Zhang           imark = l;
1334a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1335a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1336*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol));
1337a42ab0d6SHong Zhang             if (tcol) {
1338a42ab0d6SHong Zhang #else
1339a42ab0d6SHong Zhang             if ((tcol = cmap_i[cstart + cworkA[l]])) {
1340a42ab0d6SHong Zhang #endif
1341a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1342a42ab0d6SHong Zhang               if (!ijonly) {
1343*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(PetscArraycpy(mat_a,vworkA+l*bs2,bs2));
1344a42ab0d6SHong Zhang                 mat_a += bs2;
1345a42ab0d6SHong Zhang               }
134680d31651SHong Zhang               ilen++;
1347a42ab0d6SHong Zhang             }
1348a42ab0d6SHong Zhang           }
1349a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1350a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1351*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol));
1352a42ab0d6SHong Zhang             if (tcol) {
1353a42ab0d6SHong Zhang #else
1354a42ab0d6SHong Zhang             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1355a42ab0d6SHong Zhang #endif
1356a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1357a42ab0d6SHong Zhang               if (!ijonly) {
1358*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(PetscArraycpy(mat_a,vworkB+l*bs2,bs2));
1359a42ab0d6SHong Zhang                 mat_a += bs2;
1360a42ab0d6SHong Zhang               }
136180d31651SHong Zhang               ilen++;
1362a42ab0d6SHong Zhang             }
1363a42ab0d6SHong Zhang           }
1364a42ab0d6SHong Zhang         } else { /* allcolumns */
1365a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1366a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1367a42ab0d6SHong Zhang               *mat_j++ = ctmp;
1368*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscArraycpy(mat_a,vworkB+l*bs2,bs2));
1369a42ab0d6SHong Zhang               mat_a   += bs2;
137080d31651SHong Zhang               ilen++;
1371a42ab0d6SHong Zhang             } else break;
1372a42ab0d6SHong Zhang           }
1373a42ab0d6SHong Zhang           imark = l;
1374a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1375a42ab0d6SHong Zhang             *mat_j++ = cstart+cworkA[l];
1376a42ab0d6SHong Zhang             if (!ijonly) {
1377*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscArraycpy(mat_a,vworkA+l*bs2,bs2));
1378a42ab0d6SHong Zhang               mat_a += bs2;
1379a42ab0d6SHong Zhang             }
138080d31651SHong Zhang             ilen++;
1381a42ab0d6SHong Zhang           }
1382a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1383a42ab0d6SHong Zhang             *mat_j++ = bmap[cworkB[l]];
1384a42ab0d6SHong Zhang             if (!ijonly) {
1385*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscArraycpy(mat_a,vworkB+l*bs2,bs2));
1386a42ab0d6SHong Zhang               mat_a += bs2;
1387a42ab0d6SHong Zhang             }
138880d31651SHong Zhang             ilen++;
1389a42ab0d6SHong Zhang           }
1390a42ab0d6SHong Zhang         }
139180d31651SHong Zhang         imat_ilen[row] = ilen;
139217df9f7cSHong Zhang       }
139317df9f7cSHong Zhang     }
139417df9f7cSHong Zhang   }
139517df9f7cSHong Zhang 
139617df9f7cSHong Zhang   /* Now assemble the off proc rows */
13975dd0c08cSHong Zhang   if (!ijonly) {
1398*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrqs,r_waits4,MPI_STATUSES_IGNORE));
13995dd0c08cSHong Zhang   }
140017df9f7cSHong Zhang   for (tmp2=0; tmp2<nrqs; tmp2++) {
140117df9f7cSHong Zhang     sbuf1_i = sbuf1[pa[tmp2]];
140217df9f7cSHong Zhang     jmax    = sbuf1_i[0];
140317df9f7cSHong Zhang     ct1     = 2*jmax + 1;
140417df9f7cSHong Zhang     ct2     = 0;
140517df9f7cSHong Zhang     rbuf2_i = rbuf2[tmp2];
140617df9f7cSHong Zhang     rbuf3_i = rbuf3[tmp2];
14075dd0c08cSHong Zhang     if (!ijonly) rbuf4_i = rbuf4[tmp2];
140817df9f7cSHong Zhang     for (j=1; j<=jmax; j++) {
140917df9f7cSHong Zhang       is_no     = sbuf1_i[2*j-1];
141017df9f7cSHong Zhang       rmap_i    = rmap[is_no];
141117df9f7cSHong Zhang       if (!allcolumns[is_no]) cmap_i = cmap[is_no];
141217df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ*)submats[is_no]->data;
141317df9f7cSHong Zhang       imat_ilen = subc->ilen;
141417df9f7cSHong Zhang       imat_j    = subc->j;
141517df9f7cSHong Zhang       imat_i    = subc->i;
14165dd0c08cSHong Zhang       if (!ijonly) imat_a    = subc->a;
141717df9f7cSHong Zhang       max1      = sbuf1_i[2*j];
14185dd0c08cSHong Zhang       for (k=0; k<max1; k++,ct1++) { /* for each recved block row */
141917df9f7cSHong Zhang         row = sbuf1_i[ct1];
1420059f6b74SHong Zhang 
1421059f6b74SHong Zhang         if (allrows[is_no]) {
1422059f6b74SHong Zhang           row = sbuf1_i[ct1];
1423059f6b74SHong Zhang         } else {
142417df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1425*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscTableFind(rmap_i,row+1,&row));
142617df9f7cSHong Zhang           row--;
14272c71b3e2SJacob Faibussowitsch           PetscCheckFalse(row < 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
142817df9f7cSHong Zhang #else
142917df9f7cSHong Zhang           row = rmap_i[row];
143017df9f7cSHong Zhang #endif
1431059f6b74SHong Zhang         }
143217df9f7cSHong Zhang         ilen  = imat_ilen[row];
143317df9f7cSHong Zhang         mat_i = imat_i[row];
14345dd0c08cSHong Zhang         if (!ijonly) mat_a = imat_a + mat_i*bs2;
143517df9f7cSHong Zhang         mat_j = imat_j + mat_i;
143617df9f7cSHong Zhang         max2  = rbuf2_i[ct1];
143717df9f7cSHong Zhang         if (!allcolumns[is_no]) {
143817df9f7cSHong Zhang           for (l=0; l<max2; l++,ct2++) {
143917df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1440*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol));
144117df9f7cSHong Zhang #else
144217df9f7cSHong Zhang             tcol = cmap_i[rbuf3_i[ct2]];
144317df9f7cSHong Zhang #endif
144417df9f7cSHong Zhang             if (tcol) {
144517df9f7cSHong Zhang               *mat_j++ = tcol - 1;
14465dd0c08cSHong Zhang               if (!ijonly) {
1447*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(PetscArraycpy(mat_a,rbuf4_i+ct2*bs2,bs2));
14485dd0c08cSHong Zhang                 mat_a += bs2;
14495dd0c08cSHong Zhang               }
145017df9f7cSHong Zhang               ilen++;
145117df9f7cSHong Zhang             }
145217df9f7cSHong Zhang           }
145317df9f7cSHong Zhang         } else { /* allcolumns */
145417df9f7cSHong Zhang           for (l=0; l<max2; l++,ct2++) {
145517df9f7cSHong Zhang             *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */
14565dd0c08cSHong Zhang             if (!ijonly) {
1457*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscArraycpy(mat_a,rbuf4_i+ct2*bs2,bs2));
14585dd0c08cSHong Zhang               mat_a += bs2;
14595dd0c08cSHong Zhang             }
146017df9f7cSHong Zhang             ilen++;
146117df9f7cSHong Zhang           }
146217df9f7cSHong Zhang         }
146317df9f7cSHong Zhang         imat_ilen[row] = ilen;
146417df9f7cSHong Zhang       }
146517df9f7cSHong Zhang     }
146617df9f7cSHong Zhang   }
146717df9f7cSHong Zhang 
146880d31651SHong Zhang   if (!iscsorted) { /* sort column indices of the rows */
146980d31651SHong Zhang     MatScalar *work;
147080d31651SHong Zhang 
1471*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(bs2,&work));
147217df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
147317df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ*)submats[i]->data;
147480d31651SHong Zhang       imat_ilen = subc->ilen;
147517df9f7cSHong Zhang       imat_j    = subc->j;
147617df9f7cSHong Zhang       imat_i    = subc->i;
14774b6ceb0dSHong Zhang       if (!ijonly) imat_a = subc->a;
147817df9f7cSHong Zhang       if (allcolumns[i]) continue;
14794b6ceb0dSHong Zhang 
148017df9f7cSHong Zhang       jmax = nrow[i];
148117df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
148280d31651SHong Zhang         mat_i = imat_i[j];
148380d31651SHong Zhang         mat_j = imat_j + mat_i;
14844b6ceb0dSHong Zhang         ilen  = imat_ilen[j];
148580d31651SHong Zhang         if (ijonly) {
1486*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSortInt(ilen,mat_j));
148780d31651SHong Zhang         } else {
14884b6ceb0dSHong Zhang           mat_a = imat_a + mat_i*bs2;
1489*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work));
149017df9f7cSHong Zhang         }
149117df9f7cSHong Zhang       }
149217df9f7cSHong Zhang     }
1493*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(work));
149480d31651SHong Zhang   }
149517df9f7cSHong Zhang 
1496bbc89d27SHong Zhang   if (!ijonly) {
1497*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(r_waits4));
1498*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrqr,s_waits4,MPI_STATUSES_IGNORE));
1499*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(s_waits4));
1500bbc89d27SHong Zhang   }
150117df9f7cSHong Zhang 
150217df9f7cSHong Zhang   /* Restore the indices */
150317df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
1504ec3c739cSHong Zhang     if (!allrows[i]) {
1505*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISRestoreIndices(isrow[i],irow+i));
1506ec3c739cSHong Zhang     }
150717df9f7cSHong Zhang     if (!allcolumns[i]) {
1508*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISRestoreIndices(iscol[i],icol+i));
150917df9f7cSHong Zhang     }
151017df9f7cSHong Zhang   }
151117df9f7cSHong Zhang 
151217df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
1513*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY));
1514*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY));
151517df9f7cSHong Zhang   }
151617df9f7cSHong Zhang 
1517*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree5(*(PetscInt***)&irow,*(PetscInt***)&icol,nrow,ncol,issorted));
1518*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree5(row2proc,cmap,rmap,allcolumns,allrows));
151917df9f7cSHong Zhang 
1520bbc89d27SHong Zhang   if (!ijonly) {
1521175c2bc5SJunchao Zhang     if (sbuf_aa) {
1522*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(sbuf_aa[0]));
1523*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(sbuf_aa));
1524175c2bc5SJunchao Zhang     }
152517df9f7cSHong Zhang 
152617df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
1527*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(rbuf4[i]));
152817df9f7cSHong Zhang     }
1529*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rbuf4));
1530bbc89d27SHong Zhang   }
15317a868f3eSHong Zhang   c->ijonly = PETSC_FALSE; /* set back to the default */
15323a40ed3dSBarry Smith   PetscFunctionReturn(0);
1533d5b485abSSatish Balay }
1534