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