xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 
2 /*
3    Routines to compute overlapping regions of a parallel MPI matrix
4   and to find submatrices that were shared across processors.
5 */
6 #include <../src/mat/impls/baij/mpi/mpibaij.h>
7 #include <petscbt.h>
8 
9 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char **,PetscInt*,PetscInt**);
10 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt **,PetscInt**,PetscInt*);
11 extern PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
12 extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ"
16 PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov)
17 {
18   PetscErrorCode ierr;
19   PetscInt       i,N=C->cmap->N, bs=C->rmap->bs;
20   IS             *is_new;
21 
22   PetscFunctionBegin;
23   ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr);
24   /* Convert the indices into block format */
25   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr);
26   if (ov < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
27   for (i=0; i<ov; ++i) {
28     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr);
29   }
30   for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);}
31   ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr);
32   for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);}
33   ierr = PetscFree(is_new);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 /*
38   Sample message format:
39   If a processor A wants processor B to process some elements corresponding
40   to index sets is[1], is[5]
41   mesg [0] = 2   (no of index sets in the mesg)
42   -----------
43   mesg [1] = 1 => is[1]
44   mesg [2] = sizeof(is[1]);
45   -----------
46   mesg [5] = 5  => is[5]
47   mesg [6] = sizeof(is[5]);
48   -----------
49   mesg [7]
50   mesg [n]  data(is[1])
51   -----------
52   mesg[n+1]
53   mesg[m]  data(is[5])
54   -----------
55 
56   Notes:
57   nrqs - no of requests sent (or to be sent out)
58   nrqr - no of requests recieved (which have to be or which have been processed
59 */
60 #undef __FUNCT__
61 #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once"
62 PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[])
63 {
64   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
65   const PetscInt **idx,*idx_i;
66   PetscInt       *n,*w3,*w4,**data,len;
67   PetscErrorCode ierr;
68   PetscMPIInt    size,rank,tag1,tag2,*w2,*w1,nrqr;
69   PetscInt       Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr;
70   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p;
71   PetscMPIInt    *onodes1,*olengths1,*onodes2,*olengths2;
72   PetscBT        *table;
73   MPI_Comm       comm;
74   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
75   MPI_Status     *s_status,*recv_status;
76   char           *t_p;
77 
78   PetscFunctionBegin;
79   comm   = ((PetscObject)C)->comm;
80   size   = c->size;
81   rank   = c->rank;
82   Mbs    = c->Mbs;
83 
84   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
85   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
86 
87   ierr = PetscMalloc2(imax+1,const PetscInt*,&idx,imax,PetscInt,&n);CHKERRQ(ierr);
88 
89   for (i=0; i<imax; i++) {
90     ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr);
91     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
92   }
93 
94   /* evaluate communication - mesg to who,length of mesg, and buffer space
95      required. Based on this, buffers are allocated, and data copied into them*/
96   ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr);
97   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
98   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
99   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
100   for (i=0; i<imax; i++) {
101     ierr  = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
102     idx_i = idx[i];
103     len   = n[i];
104     for (j=0; j<len; j++) {
105       row  = idx_i[j];
106       if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries");
107       ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr);
108       w4[proc]++;
109     }
110     for (j=0; j<size; j++) {
111       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
112     }
113   }
114 
115   nrqs     = 0;              /* no of outgoing messages */
116   msz      = 0;              /* total mesg length (for all proc */
117   w1[rank] = 0;              /* no mesg sent to itself */
118   w3[rank] = 0;
119   for (i=0; i<size; i++) {
120     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
121   }
122   /* pa - is list of processors to communicate with */
123   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr);
124   for (i=0,j=0; i<size; i++) {
125     if (w1[i]) {pa[j] = i; j++;}
126   }
127 
128   /* Each message would have a header = 1 + 2*(no of IS) + data */
129   for (i=0; i<nrqs; i++) {
130     j      = pa[i];
131     w1[j] += w2[j] + 2*w3[j];
132     msz   += w1[j];
133   }
134 
135   /* Determine the number of messages to expect, their lengths, from from-ids */
136   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
137   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
138 
139   /* Now post the Irecvs corresponding to these messages */
140   ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr);
141 
142   /* Allocate Memory for outgoing messages */
143   ierr = PetscMalloc4(size,PetscInt*,&outdat,size,PetscInt*,&ptr,msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr);
144   ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr);
145   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
146   {
147     PetscInt *iptr = tmp,ict  = 0;
148     for (i=0; i<nrqs; i++) {
149       j         = pa[i];
150       iptr     +=  ict;
151       outdat[j] = iptr;
152       ict       = w1[j];
153     }
154   }
155 
156   /* Form the outgoing messages */
157   /*plug in the headers*/
158   for (i=0; i<nrqs; i++) {
159     j            = pa[i];
160     outdat[j][0] = 0;
161     ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
162     ptr[j]       = outdat[j] + 2*w3[j] + 1;
163   }
164 
165   /* Memory for doing local proc's work*/
166   {
167     ierr = PetscMalloc5(imax,PetscBT,&table, imax,PetscInt*,&data, imax,PetscInt,&isz,
168                         Mbs*imax,PetscInt,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,char,&t_p);CHKERRQ(ierr);
169     ierr = PetscMemzero(table,imax*sizeof(PetscBT));CHKERRQ(ierr);
170     ierr = PetscMemzero(data,imax*sizeof(PetscInt*));CHKERRQ(ierr);
171     ierr = PetscMemzero(isz,imax*sizeof(PetscInt));CHKERRQ(ierr);
172     ierr = PetscMemzero(d_p,Mbs*imax*sizeof(PetscInt));CHKERRQ(ierr);
173     ierr = PetscMemzero(t_p,(Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char));CHKERRQ(ierr);
174 
175     for (i=0; i<imax; i++) {
176       table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
177       data[i]  = d_p + (Mbs)*i;
178     }
179   }
180 
181   /* Parse the IS and update local tables and the outgoing buf with the data*/
182   {
183     PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j;
184     PetscBT  table_i;
185 
186     for (i=0; i<imax; i++) {
187       ierr    = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
188       n_i     = n[i];
189       table_i = table[i];
190       idx_i   = idx[i];
191       data_i  = data[i];
192       isz_i   = isz[i];
193       for (j=0;  j<n_i; j++) {  /* parse the indices of each IS */
194         row  = idx_i[j];
195         ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr);
196         if (proc != rank) { /* copy to the outgoing buffer */
197           ctr[proc]++;
198           *ptr[proc] = row;
199           ptr[proc]++;
200         } else { /* Update the local table */
201           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
202         }
203       }
204       /* Update the headers for the current IS */
205       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
206         if ((ctr_j = ctr[j])) {
207           outdat_j        = outdat[j];
208           k               = ++outdat_j[0];
209           outdat_j[2*k]   = ctr_j;
210           outdat_j[2*k-1] = i;
211         }
212       }
213       isz[i] = isz_i;
214     }
215   }
216 
217   /*  Now  post the sends */
218   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
219   for (i=0; i<nrqs; ++i) {
220     j    = pa[i];
221     ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr);
222   }
223 
224   /* No longer need the original indices*/
225   for (i=0; i<imax; ++i) {
226     ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr);
227   }
228   ierr = PetscFree2(idx,n);CHKERRQ(ierr);
229 
230   for (i=0; i<imax; ++i) {
231     ierr = ISDestroy(&is[i]);CHKERRQ(ierr);
232   }
233 
234   /* Do Local work*/
235   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
236 
237   /* Receive messages*/
238   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr);
239   if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);}
240 
241   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr);
242   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);}
243 
244   /* Phase 1 sends are complete - deallocate buffers */
245   ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr);
246   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
247 
248   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&xdata);CHKERRQ(ierr);
249   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt),&isz1);CHKERRQ(ierr);
250   ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
251   ierr = PetscFree(rbuf[0]);CHKERRQ(ierr);
252   ierr = PetscFree(rbuf);CHKERRQ(ierr);
253 
254   /* Send the data back*/
255   /* Do a global reduction to know the buffer space req for incoming messages*/
256   {
257     PetscMPIInt *rw1;
258 
259     ierr = PetscMalloc(size*sizeof(PetscInt),&rw1);CHKERRQ(ierr);
260     ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr);
261 
262     for (i=0; i<nrqr; ++i) {
263       proc      = recv_status[i].MPI_SOURCE;
264       if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch");
265       rw1[proc] = isz1[i];
266     }
267 
268     ierr = PetscFree(onodes1);CHKERRQ(ierr);
269     ierr = PetscFree(olengths1);CHKERRQ(ierr);
270 
271     /* Determine the number of messages to expect, their lengths, from from-ids */
272     ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr);
273     ierr = PetscFree(rw1);CHKERRQ(ierr);
274   }
275   /* Now post the Irecvs corresponding to these messages */
276   ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr);
277 
278   /*  Now  post the sends */
279   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
280   for (i=0; i<nrqr; ++i) {
281     j    = recv_status[i].MPI_SOURCE;
282     ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr);
283   }
284 
285   /* receive work done on other processors*/
286   {
287     PetscMPIInt idex;
288     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
289     PetscBT     table_i;
290     MPI_Status  *status2;
291 
292     ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr);
293     for (i=0; i<nrqs; ++i) {
294       ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr);
295       /* Process the message*/
296       rbuf2_i = rbuf2[idex];
297       ct1     = 2*rbuf2_i[0]+1;
298       jmax    = rbuf2[idex][0];
299       for (j=1; j<=jmax; j++) {
300         max     = rbuf2_i[2*j];
301         is_no   = rbuf2_i[2*j-1];
302         isz_i   = isz[is_no];
303         data_i  = data[is_no];
304         table_i = table[is_no];
305         for (k=0; k<max; k++,ct1++) {
306           row = rbuf2_i[ct1];
307           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
308         }
309         isz[is_no] = isz_i;
310       }
311     }
312     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);}
313     ierr = PetscFree(status2);CHKERRQ(ierr);
314   }
315 
316   for (i=0; i<imax; ++i) {
317     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
318   }
319 
320 
321   ierr = PetscFree(onodes2);CHKERRQ(ierr);
322   ierr = PetscFree(olengths2);CHKERRQ(ierr);
323 
324   ierr = PetscFree(pa);CHKERRQ(ierr);
325   ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr);
326   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
327   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
328   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
329   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
330   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
331   ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr);
332   ierr = PetscFree(s_status);CHKERRQ(ierr);
333   ierr = PetscFree(recv_status);CHKERRQ(ierr);
334   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
335   ierr = PetscFree(xdata);CHKERRQ(ierr);
336   ierr = PetscFree(isz1);CHKERRQ(ierr);
337   PetscFunctionReturn(0);
338 }
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local"
342 /*
343    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
344        the work on the local processor.
345 
346      Inputs:
347       C      - MAT_MPIBAIJ;
348       imax - total no of index sets processed at a time;
349       table  - an array of char - size = Mbs bits.
350 
351      Output:
352       isz    - array containing the count of the solution elements corresponding
353                to each index set;
354       data   - pointer to the solutions
355 */
356 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data)
357 {
358   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
359   Mat         A = c->A,B = c->B;
360   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
361   PetscInt    start,end,val,max,rstart,cstart,*ai,*aj;
362   PetscInt    *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
363   PetscBT     table_i;
364 
365   PetscFunctionBegin;
366   rstart = c->rstartbs;
367   cstart = c->cstartbs;
368   ai     = a->i;
369   aj     = a->j;
370   bi     = b->i;
371   bj     = b->j;
372   garray = c->garray;
373 
374 
375   for (i=0; i<imax; i++) {
376     data_i  = data[i];
377     table_i = table[i];
378     isz_i   = isz[i];
379     for (j=0,max=isz[i]; j<max; j++) {
380       row   = data_i[j] - rstart;
381       start = ai[row];
382       end   = ai[row+1];
383       for (k=start; k<end; k++) { /* Amat */
384         val = aj[k] + cstart;
385         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
386       }
387       start = bi[row];
388       end   = bi[row+1];
389       for (k=start; k<end; k++) { /* Bmat */
390         val = garray[bj[k]];
391         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
392       }
393     }
394     isz[i] = isz_i;
395   }
396   PetscFunctionReturn(0);
397 }
398 #undef __FUNCT__
399 #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive"
400 /*
401       MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages,
402          and return the output
403 
404          Input:
405            C    - the matrix
406            nrqr - no of messages being processed.
407            rbuf - an array of pointers to the recieved requests
408 
409          Output:
410            xdata - array of messages to be sent back
411            isz1  - size of each message
412 
413   For better efficiency perhaps we should malloc separately each xdata[i],
414 then if a remalloc is required we need only copy the data for that one row
415 rather than all previous rows as it is now where a single large chunck of
416 memory is used.
417 
418 */
419 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
420 {
421   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
422   Mat            A = c->A,B = c->B;
423   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
424   PetscErrorCode ierr;
425   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
426   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
427   PetscInt       val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
428   PetscInt       *rbuf_i,kmax,rbuf_0;
429   PetscBT        xtable;
430 
431   PetscFunctionBegin;
432   Mbs    = c->Mbs;
433   rstart = c->rstartbs;
434   cstart = c->cstartbs;
435   ai     = a->i;
436   aj     = a->j;
437   bi     = b->i;
438   bj     = b->j;
439   garray = c->garray;
440 
441 
442   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
443     rbuf_i  =  rbuf[i];
444     rbuf_0  =  rbuf_i[0];
445     ct     += rbuf_0;
446     for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; }
447   }
448 
449   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
450   else        max1 = 1;
451   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
452   ierr         = PetscMalloc(mem_estimate*sizeof(PetscInt),&xdata[0]);CHKERRQ(ierr);
453   ++no_malloc;
454   ierr         = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr);
455   ierr         = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr);
456 
457   ct3 = 0;
458   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
459     rbuf_i =  rbuf[i];
460     rbuf_0 =  rbuf_i[0];
461     ct1    =  2*rbuf_0+1;
462     ct2    =  ct1;
463     ct3    += ct1;
464     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
465       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
466       oct2 = ct2;
467       kmax = rbuf_i[2*j];
468       for (k=0; k<kmax; k++,ct1++) {
469         row = rbuf_i[ct1];
470         if (!PetscBTLookupSet(xtable,row)) {
471           if (!(ct3 < mem_estimate)) {
472             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
473             ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
474             ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
475             ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
476             xdata[0]     = tmp;
477             mem_estimate = new_estimate; ++no_malloc;
478             for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
479           }
480           xdata[i][ct2++] = row;
481           ct3++;
482         }
483       }
484       for (k=oct2,max2=ct2; k<max2; k++)  {
485         row   = xdata[i][k] - rstart;
486         start = ai[row];
487         end   = ai[row+1];
488         for (l=start; l<end; l++) {
489           val = aj[l] + cstart;
490           if (!PetscBTLookupSet(xtable,val)) {
491             if (!(ct3 < mem_estimate)) {
492               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
493               ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
494               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
495               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
496               xdata[0]     = tmp;
497               mem_estimate = new_estimate; ++no_malloc;
498               for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
499             }
500             xdata[i][ct2++] = val;
501             ct3++;
502           }
503         }
504         start = bi[row];
505         end   = bi[row+1];
506         for (l=start; l<end; l++) {
507           val = garray[bj[l]];
508           if (!PetscBTLookupSet(xtable,val)) {
509             if (!(ct3 < mem_estimate)) {
510               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
511               ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
512               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
513               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
514               xdata[0]     = tmp;
515               mem_estimate = new_estimate; ++no_malloc;
516               for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
517             }
518             xdata[i][ct2++] = val;
519             ct3++;
520           }
521         }
522       }
523       /* Update the header*/
524       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
525       xdata[i][2*j-1] = rbuf_i[2*j-1];
526     }
527     xdata[i][0] = rbuf_0;
528     xdata[i+1]  = xdata[i] + ct2;
529     isz1[i]     = ct2; /* size of each message */
530   }
531   ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr);
532   ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr);
533   PetscFunctionReturn(0);
534 }
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ"
538 PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
539 {
540   IS             *isrow_new,*iscol_new;
541   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
542   PetscErrorCode ierr;
543   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs;
544   PetscBool      colflag,*allcolumns,*allrows;
545 
546   PetscFunctionBegin;
547   /* Currently, unsorted column indices will result in inverted column indices in the resulting submatrices. */
548   for (i = 0; i < ismax; ++i) {
549     PetscBool sorted;
550     ierr = ISSorted(iscol[i], &sorted);CHKERRQ(ierr);
551     if (!sorted) SETERRQ1(((PetscObject)iscol[i])->comm, PETSC_ERR_SUP, "Column index set %D not sorted", i);
552   }
553   /* The compression and expansion should be avoided. Doesn't point
554      out errors, might change the indices, hence buggey */
555   ierr = PetscMalloc2(ismax+1,IS,&isrow_new,ismax+1,IS,&iscol_new);CHKERRQ(ierr);
556   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_new);CHKERRQ(ierr);
557   ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_new);CHKERRQ(ierr);
558 
559   /* Check for special case: each processor gets entire matrix columns */
560   ierr = PetscMalloc2(ismax+1,PetscBool,&allcolumns,ismax+1,PetscBool,&allrows);CHKERRQ(ierr);
561   for (i=0; i<ismax; i++) {
562     ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr);
563     ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr);
564     if (colflag && ncol == C->cmap->N) {
565       allcolumns[i] = PETSC_TRUE;
566     } else {
567       allcolumns[i] = PETSC_FALSE;
568     }
569 
570     ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr);
571     ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr);
572     if (colflag && nrow == C->rmap->N) {
573       allrows[i] = PETSC_TRUE;
574     } else {
575       allrows[i] = PETSC_FALSE;
576     }
577   }
578 
579   /* Allocate memory to hold all the submatrices */
580   if (scall != MAT_REUSE_MATRIX) {
581     ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr);
582   }
583   /* Determine the number of stages through which submatrices are done */
584   nmax          = 20*1000000 / (c->Nbs * sizeof(PetscInt));
585   if (!nmax) nmax = 1;
586   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
587 
588   /* Make sure every processor loops through the nstages */
589   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);CHKERRQ(ierr);
590   for (i=0,pos=0; i<nstages; i++) {
591     if (pos+nmax <= ismax) max_no = nmax;
592     else if (pos == ismax) max_no = 0;
593     else                   max_no = ismax-pos;
594     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr);
595     pos += max_no;
596   }
597 
598   for (i=0; i<ismax; i++) {
599     ierr = ISDestroy(&isrow_new[i]);CHKERRQ(ierr);
600     ierr = ISDestroy(&iscol_new[i]);CHKERRQ(ierr);
601   }
602   ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr);
603   ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr);
604   PetscFunctionReturn(0);
605 }
606 
607 #if defined (PETSC_USE_CTABLE)
608 #undef __FUNCT__
609 #define __FUNCT__ "PetscGetProc"
610 PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
611 {
612   PetscInt       nGlobalNd = proc_gnode[size];
613   PetscMPIInt    fproc;
614   PetscErrorCode ierr;
615 
616   PetscFunctionBegin;
617   ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr);
618   if (fproc > size) fproc = size;
619   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
620     if (row < proc_gnode[fproc]) fproc--;
621     else                         fproc++;
622   }
623   *rank = fproc;
624   PetscFunctionReturn(0);
625 }
626 #endif
627 
628 /* -------------------------------------------------------------------------*/
629 /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
630 #undef __FUNCT__
631 #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local"
632 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats)
633 {
634   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
635   Mat            A = c->A;
636   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
637   const PetscInt **irow,**icol,*irow_i;
638   PetscInt       *nrow,*ncol,*w3,*w4,start;
639   PetscErrorCode ierr;
640   PetscMPIInt    size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc;
641   PetscInt       **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row;
642   PetscInt       nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol;
643   PetscInt       **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
644   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
645   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
646   PetscInt       bs=C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
647   PetscInt       cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
648   PetscInt       *bmap = c->garray,ctmp,rstart=c->rstartbs;
649   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3;
650   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3;
651   MPI_Comm       comm;
652   PetscBool      flag;
653   PetscMPIInt    *onodes1,*olengths1;
654   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
655   /* variables below are used for the matrix numerical values - case of !ijonly */
656   MPI_Request    *r_waits4,*s_waits4;
657   MPI_Status     *r_status4,*s_status4;
658   MatScalar      **rbuf4,**sbuf_aa,*vals,*mat_a = PETSC_NULL,*sbuf_aa_i,*vworkA = PETSC_NULL,*vworkB = PETSC_NULL;
659   MatScalar      *a_a=a->a,*b_a=b->a;
660 
661 #if defined (PETSC_USE_CTABLE)
662   PetscInt       tt;
663   PetscTable     *rmap,*cmap,rmap_i,cmap_i=PETSC_NULL;
664 #else
665   PetscInt       **cmap,*cmap_i=PETSC_NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
666 #endif
667 
668   PetscFunctionBegin;
669   comm   = ((PetscObject)C)->comm;
670   tag0   = ((PetscObject)C)->tag;
671   size   = c->size;
672   rank   = c->rank;
673 
674   /* Get some new tags to keep the communication clean */
675   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
676   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
677   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
678 
679 #if defined(PETSC_USE_CTABLE)
680   ierr = PetscMalloc4(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol);CHKERRQ(ierr);
681 #else
682   ierr = PetscMalloc5(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol,Mbs+1,PetscInt,&rtable);CHKERRQ(ierr);
683   /* Create hash table for the mapping :row -> proc*/
684   for (i=0,j=0; i<size; i++) {
685     jmax = C->rmap->range[i+1]/bs;
686     for (; j<jmax; j++) {
687       rtable[j] = i;
688     }
689   }
690 #endif
691 
692   for (i=0; i<ismax; i++) {
693     if (allrows[i]) {
694       irow[i] = PETSC_NULL;
695       nrow[i] = C->rmap->N/bs;
696     } else {
697       ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
698       ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
699     }
700 
701     if (allcolumns[i]) {
702       icol[i] = PETSC_NULL;
703       ncol[i] = C->cmap->N/bs;
704     } else {
705       ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
706       ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
707     }
708   }
709 
710   /* evaluate communication - mesg to who,length of mesg,and buffer space
711      required. Based on this, buffers are allocated, and data copied into them*/
712   ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr);
713   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
714   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
715   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
716   for (i=0; i<ismax; i++) {
717     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
718     jmax   = nrow[i];
719     irow_i = irow[i];
720     for (j=0; j<jmax; j++) {
721       if (allrows[i]) {
722         row = j;
723       } else {
724         row  = irow_i[j];
725       }
726 #if defined (PETSC_USE_CTABLE)
727       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
728 #else
729       proc = rtable[row];
730 #endif
731       w4[proc]++;
732     }
733     for (j=0; j<size; j++) {
734       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
735     }
736   }
737 
738   nrqs     = 0;              /* no of outgoing messages */
739   msz      = 0;              /* total mesg length for all proc */
740   w1[rank] = 0;              /* no mesg sent to intself */
741   w3[rank] = 0;
742   for (i=0; i<size; i++) {
743     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
744   }
745   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/
746   for (i=0,j=0; i<size; i++) {
747     if (w1[i]) { pa[j] = i; j++; }
748   }
749 
750   /* Each message would have a header = 1 + 2*(no of IS) + data */
751   for (i=0; i<nrqs; i++) {
752     j     = pa[i];
753     w1[j] += w2[j] + 2* w3[j];
754     msz   += w1[j];
755   }
756 
757   /* Determine the number of messages to expect, their lengths, from from-ids */
758   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
759   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
760 
761   /* Now post the Irecvs corresponding to these messages */
762   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
763 
764   ierr = PetscFree(onodes1);CHKERRQ(ierr);
765   ierr = PetscFree(olengths1);CHKERRQ(ierr);
766 
767   /* Allocate Memory for outgoing messages */
768   ierr = PetscMalloc4(size,PetscInt*,&sbuf1,size,PetscInt*,&ptr,2*msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr);
769   ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
770   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
771   {
772     PetscInt *iptr = tmp,ict = 0;
773     for (i=0; i<nrqs; i++) {
774       j         = pa[i];
775       iptr     += ict;
776       sbuf1[j]  = iptr;
777       ict       = w1[j];
778     }
779   }
780 
781   /* Form the outgoing messages */
782   /* Initialise the header space */
783   for (i=0; i<nrqs; i++) {
784     j           = pa[i];
785     sbuf1[j][0] = 0;
786     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
787     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
788   }
789 
790   /* Parse the isrow and copy data into outbuf */
791   for (i=0; i<ismax; i++) {
792     ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
793     irow_i = irow[i];
794     jmax   = nrow[i];
795     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
796       if (allrows[i]) {
797         row = j;
798       } else {
799         row  = irow_i[j];
800       }
801 #if defined (PETSC_USE_CTABLE)
802       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
803 #else
804       proc = rtable[row];
805 #endif
806       if (proc != rank) { /* copy to the outgoing buf*/
807         ctr[proc]++;
808         *ptr[proc] = row;
809         ptr[proc]++;
810       }
811     }
812     /* Update the headers for the current IS */
813     for (j=0; j<size; j++) { /* Can Optimise this loop too */
814       if ((ctr_j = ctr[j])) {
815         sbuf1_j        = sbuf1[j];
816         k              = ++sbuf1_j[0];
817         sbuf1_j[2*k]   = ctr_j;
818         sbuf1_j[2*k-1] = i;
819       }
820     }
821   }
822 
823   /*  Now  post the sends */
824   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
825   for (i=0; i<nrqs; ++i) {
826     j = pa[i];
827     ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
828   }
829 
830   /* Post Recieves to capture the buffer size */
831   ierr     = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr);
832   ierr     = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr);
833   rbuf2[0] = tmp + msz;
834   for (i=1; i<nrqs; ++i) {
835     j        = pa[i];
836     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
837   }
838   for (i=0; i<nrqs; ++i) {
839     j    = pa[i];
840     ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
841   }
842 
843   /* Send to other procs the buf size they should allocate */
844 
845   /* Receive messages*/
846   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
847   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr);
848   ierr = PetscMalloc3(nrqr+1,PetscInt*,&sbuf2,nrqr,PetscInt,&req_size,nrqr,PetscInt,&req_source);CHKERRQ(ierr);
849   {
850     Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
851     PetscInt    *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
852 
853     for (i=0; i<nrqr; ++i) {
854       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
855       req_size[idex] = 0;
856       rbuf1_i         = rbuf1[idex];
857       start           = 2*rbuf1_i[0] + 1;
858       ierr            = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
859       ierr            = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr);
860       sbuf2_i         = sbuf2[idex];
861       for (j=start; j<end; j++) {
862         id               = rbuf1_i[j] - rstart;
863         ncols            = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
864         sbuf2_i[j]       = ncols;
865         req_size[idex] += ncols;
866       }
867       req_source[idex] = r_status1[i].MPI_SOURCE;
868       /* form the header */
869       sbuf2_i[0]   = req_size[idex];
870       for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
871       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
872     }
873   }
874   ierr = PetscFree(r_status1);CHKERRQ(ierr);
875   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
876 
877   /*  recv buffer sizes */
878   /* Receive messages*/
879   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr);
880   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr);
881   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr);
882   if (!ijonly) {
883     ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr);
884     ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr);
885   }
886 
887   for (i=0; i<nrqs; ++i) {
888     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
889     ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr);
890     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
891     if (!ijonly) {
892       ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr);
893       ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
894     }
895   }
896   ierr = PetscFree(r_status2);CHKERRQ(ierr);
897   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
898 
899   /* Wait on sends1 and sends2 */
900   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr);
901   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr);
902 
903   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
904   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
905   ierr = PetscFree(s_status1);CHKERRQ(ierr);
906   ierr = PetscFree(s_status2);CHKERRQ(ierr);
907   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
908   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
909 
910   /* Now allocate buffers for a->j, and send them off */
911   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr);
912   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
913   ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr);
914   for (i=1; i<nrqr; i++)  sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
915 
916   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr);
917   {
918      for (i=0; i<nrqr; i++) {
919       rbuf1_i   = rbuf1[i];
920       sbuf_aj_i = sbuf_aj[i];
921       ct1       = 2*rbuf1_i[0] + 1;
922       ct2       = 0;
923       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
924         kmax = rbuf1[i][2*j];
925         for (k=0; k<kmax; k++,ct1++) {
926           row    = rbuf1_i[ct1] - rstart;
927           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
928           ncols  = nzA + nzB;
929           cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
930 
931           /* load the column indices for this row into cols*/
932           cols  = sbuf_aj_i + ct2;
933           for (l=0; l<nzB; l++) {
934             if ((ctmp = bmap[cworkB[l]]) < cstart)  cols[l] = ctmp;
935             else break;
936           }
937           imark = l;
938           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
939           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
940           ct2 += ncols;
941         }
942       }
943       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
944     }
945   }
946   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr);
947   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr);
948 
949   /* Allocate buffers for a->a, and send them off */
950   if (!ijonly) {
951     ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr);
952     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
953     ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr);
954     for (i=1; i<nrqr; i++)  sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
955 
956     ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr);
957     {
958       for (i=0; i<nrqr; i++) {
959         rbuf1_i   = rbuf1[i];
960         sbuf_aa_i = sbuf_aa[i];
961         ct1       = 2*rbuf1_i[0]+1;
962         ct2       = 0;
963         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
964           kmax = rbuf1_i[2*j];
965           for (k=0; k<kmax; k++,ct1++) {
966             row    = rbuf1_i[ct1] - rstart;
967             nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
968             ncols  = nzA + nzB;
969             cworkA = a_j + a_i[row];     cworkB = b_j + b_i[row];
970             vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2;
971 
972             /* load the column values for this row into vals*/
973             vals  = sbuf_aa_i+ct2*bs2;
974             for (l=0; l<nzB; l++) {
975               if ((bmap[cworkB[l]]) < cstart) {
976                 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
977               }
978               else break;
979             }
980             imark = l;
981             for (l=0; l<nzA; l++) {
982               ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
983             }
984             for (l=imark; l<nzB; l++) {
985               ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
986             }
987             ct2 += ncols;
988           }
989         }
990         ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
991       }
992     }
993     ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr);
994     ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr);
995   }
996   ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr);
997   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
998 
999   /* Form the matrix */
1000   /* create col map: global col of C -> local col of submatrices */
1001   {
1002     const PetscInt *icol_i;
1003 #if defined (PETSC_USE_CTABLE)
1004     ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&cmap);CHKERRQ(ierr);
1005     for (i=0; i<ismax; i++) {
1006       if (!allcolumns[i]) {
1007         ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr);
1008         jmax   = ncol[i];
1009         icol_i = icol[i];
1010         cmap_i = cmap[i];
1011         for (j=0; j<jmax; j++) {
1012           ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1013         }
1014       } else {
1015         cmap[i] = PETSC_NULL;
1016       }
1017     }
1018 #else
1019     ierr    = PetscMalloc(ismax*sizeof(PetscInt*),&cmap);CHKERRQ(ierr);
1020     for (i=0; i<ismax; i++) {
1021       if (!allcolumns[i]) {
1022         ierr = PetscMalloc(c->Nbs*sizeof(PetscInt),&cmap[i]);CHKERRQ(ierr);
1023         ierr = PetscMemzero(cmap[i],c->Nbs*sizeof(PetscInt));CHKERRQ(ierr);
1024         jmax   = ncol[i];
1025         icol_i = icol[i];
1026         cmap_i = cmap[i];
1027         for (j=0; j<jmax; j++) {
1028           cmap_i[icol_i[j]] = j+1;
1029         }
1030       } else { /* allcolumns[i] */
1031         cmap[i] = PETSC_NULL;
1032       }
1033     }
1034 #endif
1035   }
1036 
1037   /* Create lens which is required for MatCreate... */
1038   for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
1039   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr);
1040   lens[0] = (PetscInt*)(lens + ismax);
1041   ierr    = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr);
1042   for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
1043 
1044   /* Update lens from local data */
1045   for (i=0; i<ismax; i++) {
1046     jmax   = nrow[i];
1047     if (!allcolumns[i]) cmap_i = cmap[i];
1048     irow_i = irow[i];
1049     lens_i = lens[i];
1050     for (j=0; j<jmax; j++) {
1051       if (allrows[i]) {
1052         row = j;
1053       } else {
1054         row  = irow_i[j];
1055       }
1056 #if defined (PETSC_USE_CTABLE)
1057       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1058 #else
1059       proc = rtable[row];
1060 #endif
1061       if (proc == rank) {
1062         /* Get indices from matA and then from matB */
1063         row    = row - rstart;
1064         nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
1065         cworkA =  a_j + a_i[row]; cworkB = b_j + b_i[row];
1066         if (!allcolumns[i]) {
1067 #if defined (PETSC_USE_CTABLE)
1068           for (k=0; k<nzA; k++) {
1069             ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1070             if (tt) { lens_i[j]++; }
1071           }
1072           for (k=0; k<nzB; k++) {
1073             ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1074             if (tt) { lens_i[j]++; }
1075           }
1076 
1077 #else
1078           for (k=0; k<nzA; k++) {
1079             if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; }
1080           }
1081           for (k=0; k<nzB; k++) {
1082             if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; }
1083           }
1084 #endif
1085         } else {/* allcolumns */
1086           lens_i[j] = nzA + nzB;
1087         }
1088       }
1089     }
1090   }
1091 #if defined (PETSC_USE_CTABLE)
1092   /* Create row map*/
1093   ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rmap);CHKERRQ(ierr);
1094   for (i=0; i<ismax; i++) {
1095     ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr);
1096   }
1097 #else
1098   /* Create row map*/
1099   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr);
1100   rmap[0] = (PetscInt*)(rmap + ismax);
1101   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr);
1102   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;}
1103 #endif
1104   for (i=0; i<ismax; i++) {
1105     irow_i = irow[i];
1106     jmax   = nrow[i];
1107 #if defined (PETSC_USE_CTABLE)
1108     rmap_i = rmap[i];
1109     for (j=0; j<jmax; j++) {
1110       if (allrows[i]) {
1111         ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1112       } else {
1113         ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1114       }
1115     }
1116 #else
1117     rmap_i = rmap[i];
1118     for (j=0; j<jmax; j++) {
1119       if (allrows[i]) {
1120         rmap_i[j] = j;
1121       } else {
1122         rmap_i[irow_i[j]] = j;
1123       }
1124     }
1125 #endif
1126   }
1127 
1128   /* Update lens from offproc data */
1129   {
1130     PetscInt    *rbuf2_i,*rbuf3_i,*sbuf1_i;
1131     PetscMPIInt ii;
1132 
1133     for (tmp2=0; tmp2<nrqs; tmp2++) {
1134       ierr    = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr);
1135       idex   = pa[ii];
1136       sbuf1_i = sbuf1[idex];
1137       jmax    = sbuf1_i[0];
1138       ct1     = 2*jmax+1;
1139       ct2     = 0;
1140       rbuf2_i = rbuf2[ii];
1141       rbuf3_i = rbuf3[ii];
1142       for (j=1; j<=jmax; j++) {
1143         is_no   = sbuf1_i[2*j-1];
1144         max1    = sbuf1_i[2*j];
1145         lens_i  = lens[is_no];
1146         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1147         rmap_i = rmap[is_no];
1148         for (k=0; k<max1; k++,ct1++) {
1149 #if defined (PETSC_USE_CTABLE)
1150           ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1151           row--;
1152           if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1153 #else
1154           row  = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
1155 #endif
1156           max2 = rbuf2_i[ct1];
1157           for (l=0; l<max2; l++,ct2++) {
1158             if (!allcolumns[is_no]) {
1159 #if defined (PETSC_USE_CTABLE)
1160               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
1161               if (tt) {
1162                 lens_i[row]++;
1163               }
1164 #else
1165               if (cmap_i[rbuf3_i[ct2]]) {
1166                 lens_i[row]++;
1167               }
1168 #endif
1169             } else { /* allcolumns */
1170               lens_i[row]++;
1171             }
1172           }
1173         }
1174       }
1175     }
1176   }
1177   ierr = PetscFree(r_status3);CHKERRQ(ierr);
1178   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
1179   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
1180   ierr = PetscFree(s_status3);CHKERRQ(ierr);
1181   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1182 
1183   /* Create the submatrices */
1184   if (scall == MAT_REUSE_MATRIX) {
1185     if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet");
1186     /*
1187         Assumes new rows are same length as the old rows, hence bug!
1188     */
1189     for (i=0; i<ismax; i++) {
1190       mat = (Mat_SeqBAIJ *)(submats[i]->data);
1191       if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1192       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
1193       if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1194       /* Initial matrix as if empty */
1195       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr);
1196       submats[i]->factortype = C->factortype;
1197     }
1198   } else {
1199     PetscInt bs_tmp;
1200     if (ijonly) {
1201       bs_tmp = 1;
1202     } else {
1203       bs_tmp = bs;
1204     }
1205     for (i=0; i<ismax; i++) {
1206       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
1207       ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr);
1208       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
1209       ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr);
1210       ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */
1211     }
1212   }
1213 
1214   /* Assemble the matrices */
1215   /* First assemble the local rows */
1216   {
1217     PetscInt  ilen_row,*imat_ilen,*imat_j,*imat_i;
1218     MatScalar *imat_a = PETSC_NULL;
1219 
1220     for (i=0; i<ismax; i++) {
1221       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1222       imat_ilen = mat->ilen;
1223       imat_j    = mat->j;
1224       imat_i    = mat->i;
1225       if (!ijonly) imat_a = mat->a;
1226       if (!allcolumns[i]) cmap_i = cmap[i];
1227       rmap_i = rmap[i];
1228       irow_i    = irow[i];
1229       jmax      = nrow[i];
1230       for (j=0; j<jmax; j++) {
1231         if (allrows[i]) {
1232           row = j;
1233         } else {
1234           row      = irow_i[j];
1235         }
1236 #if defined (PETSC_USE_CTABLE)
1237         ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1238 #else
1239         proc = rtable[row];
1240 #endif
1241         if (proc == rank) {
1242           row      = row - rstart;
1243           nzA      = a_i[row+1] - a_i[row];
1244           nzB      = b_i[row+1] - b_i[row];
1245           cworkA   = a_j + a_i[row];
1246           cworkB   = b_j + b_i[row];
1247           if (!ijonly) {
1248             vworkA = a_a + a_i[row]*bs2;
1249             vworkB = b_a + b_i[row]*bs2;
1250           }
1251 #if defined (PETSC_USE_CTABLE)
1252           ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr);
1253           row--;
1254           if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1255 #else
1256           row      = rmap_i[row + rstart];
1257 #endif
1258           mat_i    = imat_i[row];
1259           if (!ijonly) mat_a = imat_a + mat_i*bs2;
1260           mat_j    = imat_j + mat_i;
1261           ilen_row = imat_ilen[row];
1262 
1263           /* load the column indices for this row into cols*/
1264           if (!allcolumns[i]) {
1265             for (l=0; l<nzB; l++) {
1266               if ((ctmp = bmap[cworkB[l]]) < cstart) {
1267 #if defined (PETSC_USE_CTABLE)
1268                 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr);
1269                 if (tcol) {
1270 #else
1271                 if ((tcol = cmap_i[ctmp])) {
1272 #endif
1273                   *mat_j++ = tcol - 1;
1274                   ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1275                   mat_a   += bs2;
1276                   ilen_row++;
1277                 }
1278               } else break;
1279             }
1280             imark = l;
1281             for (l=0; l<nzA; l++) {
1282 #if defined (PETSC_USE_CTABLE)
1283               ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1284               if (tcol) {
1285 #else
1286               if ((tcol = cmap_i[cstart + cworkA[l]])) {
1287 #endif
1288                 *mat_j++ = tcol - 1;
1289                 if (!ijonly) {
1290                   ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1291                   mat_a += bs2;
1292                 }
1293                 ilen_row++;
1294               }
1295             }
1296             for (l=imark; l<nzB; l++) {
1297 #if defined (PETSC_USE_CTABLE)
1298               ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1299               if (tcol) {
1300 #else
1301               if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1302 #endif
1303                 *mat_j++ = tcol - 1;
1304                 if (!ijonly) {
1305                   ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1306                   mat_a   += bs2;
1307                 }
1308                 ilen_row++;
1309               }
1310             }
1311           } else { /* allcolumns */
1312             for (l=0; l<nzB; l++) {
1313               if ((ctmp = bmap[cworkB[l]]) < cstart) {
1314                 *mat_j++ = ctmp;
1315                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1316                 mat_a   += bs2;
1317                 ilen_row++;
1318               } else break;
1319             }
1320             imark = l;
1321             for (l=0; l<nzA; l++) {
1322               *mat_j++ = cstart+cworkA[l];
1323               if (!ijonly) {
1324                 ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1325                 mat_a += bs2;
1326               }
1327               ilen_row++;
1328             }
1329             for (l=imark; l<nzB; l++) {
1330               *mat_j++ = bmap[cworkB[l]];
1331               if (!ijonly) {
1332                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1333                 mat_a   += bs2;
1334               }
1335               ilen_row++;
1336             }
1337           }
1338           imat_ilen[row] = ilen_row;
1339         }
1340       }
1341     }
1342   }
1343 
1344   /*   Now assemble the off proc rows*/
1345   {
1346     PetscInt    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1347     PetscInt    *imat_j,*imat_i;
1348     MatScalar   *imat_a = PETSC_NULL,*rbuf4_i = PETSC_NULL;
1349     PetscMPIInt ii;
1350 
1351     for (tmp2=0; tmp2<nrqs; tmp2++) {
1352       if (ijonly) {
1353         ii = tmp2;
1354       } else {
1355         ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr);
1356       }
1357       idex   = pa[ii];
1358       sbuf1_i = sbuf1[idex];
1359       jmax    = sbuf1_i[0];
1360       ct1     = 2*jmax + 1;
1361       ct2     = 0;
1362       rbuf2_i = rbuf2[ii];
1363       rbuf3_i = rbuf3[ii];
1364       if (!ijonly) rbuf4_i = rbuf4[ii];
1365       for (j=1; j<=jmax; j++) {
1366         is_no     = sbuf1_i[2*j-1];
1367         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1368         rmap_i    = rmap[is_no];
1369         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1370         imat_ilen = mat->ilen;
1371         imat_j    = mat->j;
1372         imat_i    = mat->i;
1373         if (!ijonly) imat_a = mat->a;
1374         max1      = sbuf1_i[2*j];
1375         for (k=0; k<max1; k++,ct1++) {
1376           row   = sbuf1_i[ct1];
1377 #if defined (PETSC_USE_CTABLE)
1378           ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr);
1379           row--;
1380           if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1381 #else
1382           row   = rmap_i[row];
1383 #endif
1384           ilen  = imat_ilen[row];
1385           mat_i = imat_i[row];
1386           if (!ijonly) mat_a = imat_a + mat_i*bs2;
1387           mat_j = imat_j + mat_i;
1388           max2 = rbuf2_i[ct1];
1389 
1390           if (!allcolumns[is_no]) {
1391             for (l=0; l<max2; l++,ct2++) {
1392 #if defined (PETSC_USE_CTABLE)
1393               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1394               if (tcol) {
1395 #else
1396               if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1397 #endif
1398                 *mat_j++    = tcol - 1;
1399                 if (!ijonly) {
1400                   ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1401                   mat_a += bs2;
1402                 }
1403                 ilen++;
1404               }
1405             }
1406           } else { /* allcolumns */
1407             for (l=0; l<max2; l++,ct2++) {
1408               *mat_j++    = rbuf3_i[ct2];
1409               if (!ijonly) {
1410                 ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1411                 mat_a += bs2;
1412               }
1413               ilen++;
1414             }
1415           }
1416           imat_ilen[row] = ilen;
1417         }
1418       }
1419     }
1420   }
1421   if (!ijonly) {
1422     ierr = PetscFree(r_status4);CHKERRQ(ierr);
1423     ierr = PetscFree(r_waits4);CHKERRQ(ierr);
1424     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
1425     ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1426     ierr = PetscFree(s_status4);CHKERRQ(ierr);
1427   }
1428 
1429   /* Restore the indices */
1430   for (i=0; i<ismax; i++) {
1431     if (!allrows[i]) {
1432       ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1433     }
1434     if (!allcolumns[i]) {
1435       ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1436     }
1437   }
1438 
1439   /* Destroy allocated memory */
1440 #if defined(PETSC_USE_CTABLE)
1441   ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr);
1442 #else
1443   ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr);
1444 #endif
1445   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
1446   ierr = PetscFree(pa);CHKERRQ(ierr);
1447 
1448   ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
1449   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1450   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1451   for (i=0; i<nrqr; ++i) {
1452     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1453   }
1454   for (i=0; i<nrqs; ++i) {
1455     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1456   }
1457   ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr);
1458   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1459   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1460   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
1461   if (!ijonly) {
1462     for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);}
1463     ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1464     ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1465     ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
1466   }
1467 
1468 #if defined (PETSC_USE_CTABLE)
1469   for (i=0; i<ismax; i++) {
1470     ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr);
1471   }
1472 #endif
1473   ierr = PetscFree(rmap);CHKERRQ(ierr);
1474 
1475   for (i=0; i<ismax; i++) {
1476     if (!allcolumns[i]) {
1477 #if defined (PETSC_USE_CTABLE)
1478       ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr);
1479 #else
1480       ierr = PetscFree(cmap[i]);CHKERRQ(ierr);
1481 #endif
1482     }
1483   }
1484   ierr = PetscFree(cmap);CHKERRQ(ierr);
1485   ierr = PetscFree(lens);CHKERRQ(ierr);
1486 
1487   for (i=0; i<ismax; i++) {
1488     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1489     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1490   }
1491 
1492   c->ijonly = PETSC_FALSE; /* set back to the default */
1493   PetscFunctionReturn(0);
1494 }
1495 
1496