xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision ec3c739c1873fa6b648dce142cd6f72852e3f8f0)
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,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs;
524   PetscBool      colflag,*allcolumns,*allrows;
525 
526   PetscFunctionBegin;
527   /* The compression and expansion should be avoided. Doesn't point
528      out errors, might change the indices, hence buggey */
529   ierr = PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block);CHKERRQ(ierr);
530   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block);CHKERRQ(ierr);
531   ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block);CHKERRQ(ierr);
532 
533 #if 1
534   /* Check for special case: each processor gets entire matrix columns -- rm! */
535   ierr = PetscMalloc2(ismax+1,&allcolumns,ismax+1,&allrows);CHKERRQ(ierr);
536   for (i=0; i<ismax; i++) {
537     ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr);
538     ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr);
539     if (colflag && ncol == C->cmap->N) allcolumns[i] = PETSC_TRUE;
540     else allcolumns[i] = PETSC_FALSE;
541 
542     ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr);
543     ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr);
544     if (colflag && nrow == C->rmap->N) {
545       allrows[i] = PETSC_TRUE;
546       //printf("all rows \n");
547     }
548     else allrows[i] = PETSC_FALSE;
549   }
550 #endif
551 
552   /* Allocate memory to hold all the submatrices */
553   if (scall == MAT_INITIAL_MATRIX) {
554     ierr = PetscMalloc1(ismax+1,submat);CHKERRQ(ierr);
555   }
556   /* Determine the number of stages through which submatrices are done */
557   nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt));
558   if (!nmax) nmax = 1;
559   nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0);
560 
561   /* Make sure every processor loops through the nstages */
562   ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr);
563   for (i=0,pos=0; i<nstages; i++) {
564     if (pos+nmax <= ismax) max_no = nmax;
565     else if (pos == ismax) max_no = 0;
566     else                   max_no = ismax-pos;
567 
568     PetscBool isbaij;
569     ierr = PetscObjectTypeCompare((PetscObject)C,MATMPIBAIJ,&isbaij);CHKERRQ(ierr);
570     if (isbaij) {
571       ierr = MatGetSubMatrices_MPIBAIJ_local_new(C,max_no,isrow_block+pos,iscol_block+pos,scall,*submat+pos);CHKERRQ(ierr);
572     } else {
573       ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr);
574     }
575     pos += max_no;
576   }
577 
578   for (i=0; i<ismax; i++) {
579     ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr);
580     ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr);
581   }
582   ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr);
583   ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr);
584   PetscFunctionReturn(0);
585 }
586 
587 #if defined(PETSC_USE_CTABLE)
588 PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
589 {
590   PetscInt       nGlobalNd = proc_gnode[size];
591   PetscMPIInt    fproc;
592   PetscErrorCode ierr;
593 
594   PetscFunctionBegin;
595   ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr);
596   if (fproc > size) fproc = size;
597   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
598     if (row < proc_gnode[fproc]) fproc--;
599     else                         fproc++;
600   }
601   *rank = fproc;
602   PetscFunctionReturn(0);
603 }
604 #endif
605 
606 /* -------------------------------------------------------------------------*/
607 /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
608 
609 PetscErrorCode MatDestroy_MPIBAIJ_MatGetSubmatrices(Mat C)
610 {
611   PetscErrorCode ierr;
612   Mat_SeqBAIJ    *c = (Mat_SeqBAIJ*)C->data;
613   Mat_SubMat     *submatj = c->submatis1;
614   PetscInt       i;
615 
616   PetscFunctionBegin;
617   if (!submatj->id) { /* delete data that are linked only to submats[id=0] */
618     ierr = PetscFree4(submatj->sbuf1,submatj->ptr,submatj->tmp,submatj->ctr);CHKERRQ(ierr);
619 
620     for (i=0; i<submatj->nrqr; ++i) {
621       ierr = PetscFree(submatj->sbuf2[i]);CHKERRQ(ierr);
622     }
623     ierr = PetscFree3(submatj->sbuf2,submatj->req_size,submatj->req_source1);CHKERRQ(ierr);
624 
625     if (submatj->rbuf1) {
626       ierr = PetscFree(submatj->rbuf1[0]);CHKERRQ(ierr);
627       ierr = PetscFree(submatj->rbuf1);CHKERRQ(ierr);
628     }
629 
630     for (i=0; i<submatj->nrqs; ++i) {
631       ierr = PetscFree(submatj->rbuf3[i]);CHKERRQ(ierr);
632     }
633     ierr = PetscFree3(submatj->req_source2,submatj->rbuf2,submatj->rbuf3);CHKERRQ(ierr);
634     ierr = PetscFree(submatj->pa);CHKERRQ(ierr);
635   }
636 
637 #if defined(PETSC_USE_CTABLE)
638   ierr = PetscTableDestroy((PetscTable*)&submatj->rmap);CHKERRQ(ierr);
639   if (submatj->cmap_loc) {ierr = PetscFree(submatj->cmap_loc);CHKERRQ(ierr);}
640   ierr = PetscFree(submatj->rmap_loc);CHKERRQ(ierr);
641 #else
642   ierr = PetscFree(submatj->rmap);CHKERRQ(ierr);
643 #endif
644 
645   if (!submatj->allcolumns) {
646 #if defined(PETSC_USE_CTABLE)
647     ierr = PetscTableDestroy((PetscTable*)&submatj->cmap);CHKERRQ(ierr);
648 #else
649     ierr = PetscFree(submatj->cmap);CHKERRQ(ierr);
650 #endif
651   }
652   ierr = submatj->destroy(C);CHKERRQ(ierr);
653   ierr = PetscFree(submatj->row2proc);CHKERRQ(ierr);
654 
655   ierr = PetscFree(submatj);CHKERRQ(ierr);
656   PetscFunctionReturn(0);
657 }
658 
659 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local_new(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
660 {
661   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
662   Mat            A  = c->A;
663   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc;
664   const PetscInt **icol,**irow;
665   PetscInt       *nrow,*ncol,start;
666   PetscErrorCode ierr;
667   PetscMPIInt    rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr;
668   PetscInt       **sbuf1,**sbuf2,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1;
669   PetscInt       nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol;
670   PetscInt       **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2;
671   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
672 #if defined(PETSC_USE_CTABLE)
673   PetscTable     *cmap,cmap_i=NULL,*rmap,rmap_i;
674 #else
675   PetscInt       **cmap,*cmap_i=NULL,**rmap,*rmap_i;
676 #endif
677   const PetscInt *irow_i;
678   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
679   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
680   MPI_Request    *r_waits4,*s_waits3,*s_waits4;
681   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
682   MPI_Status     *r_status3,*r_status4,*s_status4;
683   MPI_Comm       comm;
684   PetscScalar    **rbuf4,*rbuf4_i,**sbuf_aa,*vals,*mat_a,*imat_a,*sbuf_aa_i;
685   PetscMPIInt    *onodes1,*olengths1,end;
686   PetscInt       **row2proc,*row2proc_i,*imat_ilen,*imat_j,*imat_i;
687   Mat_SubMat     **smats,*smat_i;
688   PetscBool      *issorted,colflag,iscsorted=PETSC_TRUE;
689   PetscInt       *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen;
690 
691   PetscInt       bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs;
692   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
693   PetscInt       nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB;
694   PetscScalar    *vworkA,*vworkB,*a_a = a->a,*b_a = b->a;
695   PetscInt       cstart = c->cstartbs,*bmap = c->garray;
696   PetscBool      *allrows,*allcolumns;
697 
698   PetscFunctionBegin;
699   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
700   size = c->size;
701   rank = c->rank;
702 
703   ierr = PetscMalloc6(ismax,&smats,ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax,&allcolumns,ismax,&allrows);CHKERRQ(ierr);
704   ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr);
705 
706   for (i=0; i<ismax; i++) {
707     ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr);
708     if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */
709     ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr);
710 
711     /* Check for special case: allcolumns */
712     ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr);
713     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
714 
715     if (colflag && ncol[i] == c->Nbs) {
716       allcolumns[i] = PETSC_TRUE;
717       icol[i]       = NULL;
718       /* printf("[%d] allcolumns[%d] true\n",rank,i); */
719     } else {
720       allcolumns[i] = PETSC_FALSE;
721       ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
722     }
723 
724     /* Check for special case: allrows */
725     ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr);
726     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
727     if (colflag && nrow[i] == c->Mbs) {
728       allrows[i] = PETSC_TRUE;
729       irow[i]    = NULL;
730       /* printf("[%d] allrows[%d]\n",rank,i); */
731     } else {
732       allrows[i] = PETSC_FALSE;
733       ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
734     }
735   }
736 
737   if (scall == MAT_REUSE_MATRIX) {
738     /* Assumes new rows are same length as the old rows */
739     for (i=0; i<ismax; i++) {
740       subc = (Mat_SeqBAIJ*)(submats[i]->data);
741       if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
742 
743       /* Initial matrix as if empty */
744       ierr = PetscMemzero(subc->ilen,subc->mbs*sizeof(PetscInt));CHKERRQ(ierr);
745 
746       /* Initial matrix as if empty */
747       submats[i]->factortype = C->factortype;
748 
749       smat_i   = subc->submatis1;
750       smats[i] = smat_i;
751 
752       nrqs        = smat_i->nrqs;
753       nrqr        = smat_i->nrqr;
754       rbuf1       = smat_i->rbuf1;
755       rbuf2       = smat_i->rbuf2;
756       rbuf3       = smat_i->rbuf3;
757       req_source2 = smat_i->req_source2;
758 
759       sbuf1     = smat_i->sbuf1;
760       sbuf2     = smat_i->sbuf2;
761       ptr       = smat_i->ptr;
762       tmp       = smat_i->tmp;
763       ctr       = smat_i->ctr;
764 
765       pa          = smat_i->pa;
766       req_size    = smat_i->req_size;
767       req_source1 = smat_i->req_source1;
768 
769       allcolumns[i] = smat_i->allcolumns;
770       allrows[i]    = smat_i->allrows;
771       row2proc[i]   = smat_i->row2proc;
772       rmap[i]       = smat_i->rmap;
773       cmap[i]       = smat_i->cmap;
774     }
775   } else { /* scall == MAT_INITIAL_MATRIX */
776     /* Get some new tags to keep the communication clean */
777     ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
778     ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
779 
780     /* evaluate communication - mesg to who, length of mesg, and buffer space
781      required. Based on this, buffers are allocated, and data copied into them*/
782     ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr);   /* mesg size, initialize work vectors */
783 
784     for (i=0; i<ismax; i++) {
785       jmax   = nrow[i];
786       irow_i = irow[i];
787 
788       ierr   = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr);
789       row2proc[i] = row2proc_i;
790 
791       if (issorted[i]) proc = 0;
792       for (j=0; j<jmax; j++) {
793         if (!issorted[i]) proc = 0;
794         if (allrows[i]) row = j;
795         else row = irow_i[j];
796 
797         while (row >= c->rangebs[proc+1]) proc++;
798         w4[proc]++;
799         row2proc_i[j] = proc; /* map row index to proc */
800       }
801       for (j=0; j<size; j++) {
802         if (w4[j]) { w1[j] += w4[j];  w3[j]++; w4[j] = 0;}
803       }
804     }
805 
806     nrqs     = 0;              /* no of outgoing messages */
807     msz      = 0;              /* total mesg length (for all procs) */
808     w1[rank] = 0;              /* no mesg sent to self */
809     w3[rank] = 0;
810     for (i=0; i<size; i++) {
811       if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
812     }
813     ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/
814     for (i=0,j=0; i<size; i++) {
815       if (w1[i]) { pa[j] = i; j++; }
816     }
817 
818     /* Each message would have a header = 1 + 2*(no of IS) + data */
819     for (i=0; i<nrqs; i++) {
820       j      = pa[i];
821       w1[j] += w2[j] + 2* w3[j];
822       msz   += w1[j];
823     }
824     ierr = PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);CHKERRQ(ierr);
825 
826     /* Determine the number of messages to expect, their lengths, from from-ids */
827     ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
828     ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
829 
830     /* Now post the Irecvs corresponding to these messages */
831     tag0 = ((PetscObject)C)->tag;
832     ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
833     /* printf("[%d] nrqs %d, nrqr %d\n",rank,nrqs,nrqr); */
834 
835     ierr = PetscFree(onodes1);CHKERRQ(ierr);
836     ierr = PetscFree(olengths1);CHKERRQ(ierr);
837 
838     /* Allocate Memory for outgoing messages */
839     ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr);
840     ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
841     ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
842 
843     {
844       PetscInt *iptr = tmp;
845       k    = 0;
846       for (i=0; i<nrqs; i++) {
847         j        = pa[i];
848         iptr    += k;
849         sbuf1[j] = iptr;
850         k        = w1[j];
851       }
852     }
853 
854     /* Form the outgoing messages. Initialize the header space */
855     for (i=0; i<nrqs; i++) {
856       j           = pa[i];
857       sbuf1[j][0] = 0;
858       ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
859       ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
860     }
861 
862     /* Parse the isrow and copy data into outbuf */
863     for (i=0; i<ismax; i++) {
864       row2proc_i = row2proc[i];
865       ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
866       irow_i = irow[i];
867       jmax   = nrow[i];
868       for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
869         proc = row2proc_i[j];
870         if (allrows[i]) row = j;
871         else row = irow_i[j];
872 
873         if (proc != rank) { /* copy to the outgoing buf*/
874           ctr[proc]++;
875           *ptr[proc] = row;
876           ptr[proc]++;
877         }
878       }
879       /* Update the headers for the current IS */
880       for (j=0; j<size; j++) { /* Can Optimise this loop too */
881         if ((ctr_j = ctr[j])) {
882           sbuf1_j        = sbuf1[j];
883           k              = ++sbuf1_j[0];
884           sbuf1_j[2*k]   = ctr_j;
885           sbuf1_j[2*k-1] = i;
886         }
887       }
888     }
889 
890     /*  Now  post the sends */
891     ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr);
892     for (i=0; i<nrqs; ++i) {
893       j    = pa[i];
894       ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
895     }
896 
897     /* Post Receives to capture the buffer size */
898     ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr);
899     ierr = PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);CHKERRQ(ierr);
900     rbuf2[0] = tmp + msz;
901     for (i=1; i<nrqs; ++i) {
902       rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
903     }
904     for (i=0; i<nrqs; ++i) {
905       j    = pa[i];
906       ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRQ(ierr);
907     }
908 
909     /* Send to other procs the buf size they should allocate */
910     /* Receive messages*/
911     ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr);
912     ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr);
913     ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr);
914     {
915       PetscInt   *sAi = a->i,*sBi = b->i,id;
916       PetscInt   *sbuf2_i;
917 
918       ierr = MPI_Waitall(nrqr,r_waits1,r_status1);CHKERRQ(ierr);
919       for (i=0; i<nrqr; ++i) {
920         req_size[i] = 0;
921         rbuf1_i        = rbuf1[i];
922         start          = 2*rbuf1_i[0] + 1;
923         ierr           = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
924         ierr           = PetscMalloc1(end+1,&sbuf2[i]);CHKERRQ(ierr);
925         sbuf2_i        = sbuf2[i];
926         for (j=start; j<end; j++) {
927           id              = rbuf1_i[j] - rstart;
928           ncols           = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
929           sbuf2_i[j]      = ncols;
930           req_size[i] += ncols;
931         }
932         req_source1[i] = r_status1[i].MPI_SOURCE;
933         /* form the header */
934         sbuf2_i[0] = req_size[i];
935         for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
936 
937         ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRQ(ierr);
938       }
939     }
940     ierr = PetscFree(r_status1);CHKERRQ(ierr);
941     ierr = PetscFree(r_waits1);CHKERRQ(ierr);
942     ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
943 
944     /* Receive messages*/
945     ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr);
946     ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr);
947 
948     ierr = MPI_Waitall(nrqs,r_waits2,r_status2);CHKERRQ(ierr);
949     for (i=0; i<nrqs; ++i) {
950       ierr = PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);CHKERRQ(ierr);
951       req_source2[i] = r_status2[i].MPI_SOURCE;
952       ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRQ(ierr);
953     }
954     ierr = PetscFree(r_status2);CHKERRQ(ierr);
955     ierr = PetscFree(r_waits2);CHKERRQ(ierr);
956 
957     /* Wait on sends1 and sends2 */
958     ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr);
959     ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr);
960 
961     if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
962     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
963     ierr = PetscFree(s_status1);CHKERRQ(ierr);
964     ierr = PetscFree(s_status2);CHKERRQ(ierr);
965     ierr = PetscFree(s_waits1);CHKERRQ(ierr);
966     ierr = PetscFree(s_waits2);CHKERRQ(ierr);
967 
968     /* Now allocate sending buffers for a->j, and send them off */
969     ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr);
970     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
971     ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr);
972     for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
973 
974     ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr);
975     {
976 
977       for (i=0; i<nrqr; i++) {
978         rbuf1_i   = rbuf1[i];
979         sbuf_aj_i = sbuf_aj[i];
980         ct1       = 2*rbuf1_i[0] + 1;
981         ct2       = 0;
982         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
983           kmax = rbuf1[i][2*j];
984           for (k=0; k<kmax; k++,ct1++) {
985             row    = rbuf1_i[ct1] - rstart;
986             nzA    = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row];
987             ncols  = nzA + nzB;
988             cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
989 
990             /* load the column indices for this row into cols */
991             cols = sbuf_aj_i + ct2;
992             for (l=0; l<nzB; l++) {
993               if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
994               else break;
995             }
996             imark = l;
997             for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];}
998             for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
999             ct2 += ncols;
1000           }
1001         }
1002         ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRQ(ierr);
1003       }
1004     }
1005     ierr = PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);CHKERRQ(ierr);
1006 
1007     /* create col map: global col of C -> local col of submatrices */
1008     const PetscInt *icol_i;
1009 #if defined(PETSC_USE_CTABLE)
1010     for (i=0; i<ismax; i++) {
1011       if (!allcolumns[i]) {
1012         ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr);
1013 
1014         jmax   = ncol[i];
1015         icol_i = icol[i];
1016         cmap_i = cmap[i];
1017         for (j=0; j<jmax; j++) {
1018           ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1019         }
1020       } else cmap[i] = NULL;
1021     }
1022 #else
1023     for (i=0; i<ismax; i++) {
1024       if (!allcolumns[i]) {
1025         ierr   = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr);
1026         jmax   = ncol[i];
1027         icol_i = icol[i];
1028         cmap_i = cmap[i];
1029         for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1;
1030       } else cmap[i] = NULL;
1031     }
1032 #endif
1033 
1034     /* Create lens which is required for MatCreate... */
1035     for (i=0,j=0; i<ismax; i++) j += nrow[i];
1036     ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr);
1037 
1038     if (ismax) {
1039       ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr);
1040     }
1041     for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];
1042 
1043     /* Update lens from local data */
1044     for (i=0; i<ismax; i++) {
1045       row2proc_i = row2proc[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]) row = j;
1052         else row = irow_i[j]; /* global blocked row of C */
1053 
1054         proc = row2proc_i[j];
1055         if (proc == rank) {
1056           /* Get indices from matA and then from matB */
1057 #if defined(PETSC_USE_CTABLE)
1058           PetscInt   tt;
1059 #endif
1060           row    = row - rstart;
1061           nzA    = a_i[row+1] - a_i[row];
1062           nzB    = b_i[row+1] - b_i[row];
1063           cworkA =  a_j + a_i[row];
1064           cworkB = b_j + b_i[row];
1065 
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 
1092     /* Create row map: global row of C -> local row of submatrices */
1093 #if defined(PETSC_USE_CTABLE)
1094     for (i=0; i<ismax; i++) {
1095       ierr   = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr);
1096       irow_i = irow[i];
1097       jmax   = nrow[i];
1098       for (j=0; j<jmax; j++) {
1099         if (allrows[i]) {
1100           ierr = PetscTableAdd(rmap[i],j+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1101         } else {
1102           ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1103         }
1104       }
1105     }
1106 #else
1107     for (i=0; i<ismax; i++) {
1108       ierr   = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr);
1109       rmap_i = rmap[i];
1110       irow_i = irow[i];
1111       jmax   = nrow[i];
1112       for (j=0; j<jmax; j++) {
1113         if (allrows[i]) rmap_i[j] = j;
1114         else rmap_i[irow_i[j]] = j;
1115       }
1116     }
1117 #endif
1118 
1119     /* Update lens from offproc data */
1120     {
1121       PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i;
1122 
1123       ierr    = MPI_Waitall(nrqs,r_waits3,r_status3);CHKERRQ(ierr);
1124       for (tmp2=0; tmp2<nrqs; tmp2++) {
1125         sbuf1_i = sbuf1[pa[tmp2]];
1126         jmax    = sbuf1_i[0];
1127         ct1     = 2*jmax+1;
1128         ct2     = 0;
1129         rbuf2_i = rbuf2[tmp2];
1130         rbuf3_i = rbuf3[tmp2];
1131         for (j=1; j<=jmax; j++) {
1132           is_no  = sbuf1_i[2*j-1];
1133           max1   = sbuf1_i[2*j];
1134           lens_i = lens[is_no];
1135           if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1136           rmap_i = rmap[is_no];
1137           for (k=0; k<max1; k++,ct1++) {
1138 #if defined(PETSC_USE_CTABLE)
1139             ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1140             row--;
1141             if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1142 #else
1143             row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
1144 #endif
1145             max2 = rbuf2_i[ct1];
1146             for (l=0; l<max2; l++,ct2++) {
1147               if (!allcolumns[is_no]) {
1148 #if defined(PETSC_USE_CTABLE)
1149                 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1150 #else
1151                 tcol = cmap_i[rbuf3_i[ct2]];
1152 #endif
1153                 if (tcol) lens_i[row]++;
1154               } else { /* allcolumns */
1155                 lens_i[row]++; /* lens_i[row] += max2 ? */
1156               }
1157             }
1158           }
1159         }
1160       }
1161     }
1162     ierr = PetscFree(r_waits3);CHKERRQ(ierr);
1163     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
1164     ierr = PetscFree2(r_status3,s_status3);CHKERRQ(ierr);
1165     ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1166 
1167     /* Create the submatrices */
1168     for (i=0; i<ismax; i++) {
1169       PetscInt bs_tmp;
1170       if (ijonly) bs_tmp = 1;
1171       else        bs_tmp = bs;
1172 
1173       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
1174       ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1175 
1176       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
1177       ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr);
1178       ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */
1179 
1180       /* create struct Mat_SubMat and attached it to submat */
1181       ierr = PetscNew(&smat_i);CHKERRQ(ierr);
1182       subc = (Mat_SeqBAIJ*)submats[i]->data;
1183       subc->submatis1 = smat_i;
1184       smats[i]        = smat_i;
1185 
1186       smat_i->destroy          = submats[i]->ops->destroy;
1187       submats[i]->ops->destroy = MatDestroy_MPIBAIJ_MatGetSubmatrices;
1188       submats[i]->factortype   = C->factortype;
1189 
1190       smat_i->id          = i;
1191       smat_i->nrqs        = nrqs;
1192       smat_i->nrqr        = nrqr;
1193       smat_i->rbuf1       = rbuf1;
1194       smat_i->rbuf2       = rbuf2;
1195       smat_i->rbuf3       = rbuf3;
1196       smat_i->sbuf2       = sbuf2;
1197       smat_i->req_source2 = req_source2;
1198 
1199       smat_i->sbuf1       = sbuf1;
1200       smat_i->ptr         = ptr;
1201       smat_i->tmp         = tmp;
1202       smat_i->ctr         = ctr;
1203 
1204       smat_i->pa           = pa;
1205       smat_i->req_size     = req_size;
1206       smat_i->req_source1  = req_source1;
1207 
1208       smat_i->allcolumns  = allcolumns[i];
1209       smat_i->allrows     = allrows[i];
1210       smat_i->singleis    = PETSC_FALSE;
1211       smat_i->row2proc    = row2proc[i];
1212       smat_i->rmap        = rmap[i];
1213       smat_i->cmap        = cmap[i];
1214     }
1215 
1216     if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);}
1217     ierr = PetscFree(lens);CHKERRQ(ierr);
1218     ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1219     ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
1220 
1221   } /* endof scall == MAT_INITIAL_MATRIX */
1222 
1223   /* Post recv matrix values */
1224   ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr);
1225   ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr);
1226   ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr);
1227   ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr);
1228   ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr);
1229   for (i=0; i<nrqs; ++i) {
1230     ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr);
1231     ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRQ(ierr);
1232   }
1233 
1234   /* Allocate sending buffers for a->a, and send them off */
1235   ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr);
1236   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1237 
1238   ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr);
1239   for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
1240 
1241   ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr);
1242 
1243   for (i=0; i<nrqr; i++) {
1244     rbuf1_i   = rbuf1[i];
1245     sbuf_aa_i = sbuf_aa[i];
1246     ct1       = 2*rbuf1_i[0]+1;
1247     ct2       = 0;
1248     for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
1249       kmax = rbuf1_i[2*j];
1250       for (k=0; k<kmax; k++,ct1++) {
1251         row    = rbuf1_i[ct1] - rstart;
1252         nzA    = a_i[row+1] - a_i[row];
1253         nzB    = b_i[row+1] - b_i[row];
1254         ncols  = nzA + nzB;
1255         cworkB = b_j + b_i[row];
1256         vworkA = a_a + a_i[row]*bs2;
1257         vworkB = b_a + b_i[row]*bs2;
1258 
1259         /* load the column values for this row into vals*/
1260         vals = sbuf_aa_i+ct2*bs2;
1261         for (l=0; l<nzB; l++) {
1262           if ((bmap[cworkB[l]]) < cstart) {
1263             ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1264           } else break;
1265         }
1266         imark = l;
1267         for (l=0; l<nzA; l++) {
1268           ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);//error in proc[1]???
1269         }
1270         for (l=imark; l<nzB; l++) {
1271           ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1272         }
1273 
1274         ct2 += ncols;
1275       }
1276     }
1277     ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRQ(ierr);
1278   }
1279 
1280   if (!ismax) {
1281     ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr);
1282     ierr = PetscFree(rbuf1);CHKERRQ(ierr);
1283   }
1284 
1285   /* Assemble the matrices */
1286   /* First assemble the local rows */
1287   for (i=0; i<ismax; i++) {
1288     row2proc_i = row2proc[i];
1289     subc      = (Mat_SeqBAIJ*)submats[i]->data;
1290     imat_ilen = subc->ilen;
1291     imat_j    = subc->j;
1292     imat_i    = subc->i;
1293     imat_a    = subc->a;
1294 
1295     if (!allcolumns[i]) cmap_i = cmap[i];
1296     rmap_i = rmap[i];
1297     irow_i = irow[i];
1298     jmax   = nrow[i];
1299     for (j=0; j<jmax; j++) {
1300       if (allrows[i]) row = j;
1301       else row  = irow_i[j];
1302       proc = row2proc_i[j];
1303 
1304       if (proc == rank) {
1305 
1306         row    = row - rstart;
1307         nzA    = a_i[row+1] - a_i[row];
1308         nzB    = b_i[row+1] - b_i[row];
1309         cworkA = a_j + a_i[row];
1310         cworkB = b_j + b_i[row];
1311         if (!ijonly) {
1312           vworkA = a_a + a_i[row]*bs2;
1313           vworkB = b_a + b_i[row]*bs2;
1314         }
1315 #if defined(PETSC_USE_CTABLE)
1316         ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr);
1317         row--;
1318         if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1319 #else
1320         row = rmap_i[row + rstart];
1321 #endif
1322         mat_i = imat_i[row];
1323         if (!ijonly) mat_a = imat_a + mat_i*bs2;
1324         mat_j    = imat_j + mat_i;
1325         ilen = imat_ilen[row];
1326 
1327         /* load the column indices for this row into cols*/
1328         if (!allcolumns[i]) {
1329           for (l=0; l<nzB; l++) {
1330             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1331 #if defined(PETSC_USE_CTABLE)
1332               ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr);
1333               if (tcol) {
1334 #else
1335               if ((tcol = cmap_i[ctmp])) {
1336 #endif
1337                 *mat_j++ = tcol - 1;
1338                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1339                 mat_a   += bs2;
1340                 ilen++;
1341               }
1342             } else break;
1343           }
1344           imark = l;
1345           for (l=0; l<nzA; l++) {
1346 #if defined(PETSC_USE_CTABLE)
1347             ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1348             if (tcol) {
1349 #else
1350             if ((tcol = cmap_i[cstart + cworkA[l]])) {
1351 #endif
1352               *mat_j++ = tcol - 1;
1353               if (!ijonly) {
1354                 ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1355                 mat_a += bs2;
1356               }
1357               ilen++;
1358             }
1359           }
1360           for (l=imark; l<nzB; l++) {
1361 #if defined(PETSC_USE_CTABLE)
1362             ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1363             if (tcol) {
1364 #else
1365             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1366 #endif
1367               *mat_j++ = tcol - 1;
1368               if (!ijonly) {
1369                 ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1370                 mat_a += bs2;
1371               }
1372               ilen++;
1373             }
1374           }
1375         } else { /* allcolumns */
1376           for (l=0; l<nzB; l++) {
1377             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1378               *mat_j++ = ctmp;
1379               ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1380               mat_a   += bs2;
1381               ilen++;
1382             } else break;
1383           }
1384           imark = l;
1385           for (l=0; l<nzA; l++) {
1386             *mat_j++ = cstart+cworkA[l];
1387             if (!ijonly) {
1388               ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1389               mat_a += bs2;
1390             }
1391             ilen++;
1392           }
1393           for (l=imark; l<nzB; l++) {
1394             *mat_j++ = bmap[cworkB[l]];
1395             if (!ijonly) {
1396               ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1397               mat_a += bs2;
1398             }
1399             ilen++;
1400           }
1401         }
1402         imat_ilen[row] = ilen;
1403       }
1404     }
1405   }
1406 
1407   /* Now assemble the off proc rows */
1408   if (!ijonly) {
1409     ierr = MPI_Waitall(nrqs,r_waits4,r_status4);CHKERRQ(ierr);
1410   }
1411   for (tmp2=0; tmp2<nrqs; tmp2++) {
1412     sbuf1_i = sbuf1[pa[tmp2]];
1413     jmax    = sbuf1_i[0];
1414     ct1     = 2*jmax + 1;
1415     ct2     = 0;
1416     rbuf2_i = rbuf2[tmp2];
1417     rbuf3_i = rbuf3[tmp2];
1418     if (!ijonly) rbuf4_i = rbuf4[tmp2];
1419     for (j=1; j<=jmax; j++) {
1420       is_no     = sbuf1_i[2*j-1];
1421       rmap_i    = rmap[is_no];
1422       if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1423       subc      = (Mat_SeqBAIJ*)submats[is_no]->data;
1424       imat_ilen = subc->ilen;
1425       imat_j    = subc->j;
1426       imat_i    = subc->i;
1427       if (!ijonly) imat_a    = subc->a;
1428       max1      = sbuf1_i[2*j];
1429       for (k=0; k<max1; k++,ct1++) { /* for each recved block row */
1430         row = sbuf1_i[ct1];
1431 #if defined(PETSC_USE_CTABLE)
1432         ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr);
1433         row--;
1434         if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1435 #else
1436         row = rmap_i[row];
1437 #endif
1438         ilen  = imat_ilen[row];
1439         mat_i = imat_i[row];
1440         if (!ijonly) mat_a = imat_a + mat_i*bs2;
1441         mat_j = imat_j + mat_i;
1442         max2  = rbuf2_i[ct1];
1443         if (!allcolumns[is_no]) {
1444           for (l=0; l<max2; l++,ct2++) {
1445 #if defined(PETSC_USE_CTABLE)
1446             ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1447 #else
1448             tcol = cmap_i[rbuf3_i[ct2]];
1449 #endif
1450             if (tcol) {
1451               *mat_j++ = tcol - 1;
1452               if (!ijonly) {
1453                 ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1454                 mat_a += bs2;
1455               }
1456               //*mat_a++ = rbuf4_i[ct2];
1457               ilen++;
1458             }
1459           }
1460         } else { /* allcolumns */
1461           for (l=0; l<max2; l++,ct2++) {
1462             *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */
1463             //*mat_a++ = rbuf4_i[ct2];
1464             if (!ijonly) {
1465               ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1466               mat_a += bs2;
1467             }
1468             ilen++;
1469           }
1470         }
1471         imat_ilen[row] = ilen;
1472       }
1473     }
1474   }
1475 
1476   if (!iscsorted) { /* sort column indices of the rows */
1477     MatScalar *work;
1478 
1479     ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr);
1480     for (i=0; i<ismax; i++) {
1481       subc      = (Mat_SeqBAIJ*)submats[i]->data;
1482       imat_ilen = subc->ilen;
1483       imat_j    = subc->j;
1484       imat_i    = subc->i;
1485       if (ijonly) imat_a = subc->a;
1486 
1487       if (allcolumns[i]) continue;
1488       jmax = nrow[i];
1489       for (j=0; j<jmax; j++) {
1490         ilen  = imat_ilen[j];
1491         mat_i = imat_i[j];
1492         mat_j = imat_j + mat_i;
1493         if (ijonly) {
1494           ierr = PetscSortInt(ilen,mat_j);CHKERRQ(ierr);
1495         } else {
1496           mat_a = imat_a + mat_i;
1497           ierr  = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);
1498         }
1499       }
1500     }
1501     ierr = PetscFree(work);CHKERRQ(ierr);
1502   }
1503 
1504   ierr = PetscFree(r_status4);CHKERRQ(ierr);
1505   ierr = PetscFree(r_waits4);CHKERRQ(ierr);
1506   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
1507   ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1508   ierr = PetscFree(s_status4);CHKERRQ(ierr);
1509 
1510   /* Restore the indices */
1511   for (i=0; i<ismax; i++) {
1512     if (!allrows[i]) {
1513       ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1514     }
1515     if (!allcolumns[i]) {
1516       ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1517     }
1518   }
1519 
1520   for (i=0; i<ismax; i++) {
1521     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1522     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1523   }
1524 
1525   /* Destroy allocated memory */
1526   if (!ismax) {
1527     ierr = PetscFree(pa);CHKERRQ(ierr);
1528 
1529     ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
1530     for (i=0; i<nrqr; ++i) {
1531       ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1532     }
1533     for (i=0; i<nrqs; ++i) {
1534       ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1535     }
1536 
1537     ierr = PetscFree3(sbuf2,req_size,req_source1);CHKERRQ(ierr);
1538     ierr = PetscFree3(req_source2,rbuf2,rbuf3);CHKERRQ(ierr);
1539   }
1540 
1541   ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1542   ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
1543   ierr = PetscFree5(irow,icol,nrow,ncol,issorted);CHKERRQ(ierr);
1544 
1545   for (i=0; i<nrqs; ++i) {
1546     ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);
1547   }
1548   ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1549 
1550   ierr = PetscFree6(smats,row2proc,cmap,rmap,allcolumns,allrows);CHKERRQ(ierr);
1551   PetscFunctionReturn(0);
1552 }
1553 
1554 //------------ endof new -------------
1555 
1556 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats)
1557 {
1558   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
1559   Mat            A  = c->A;
1560   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
1561   const PetscInt **irow,**icol,*irow_i;
1562   PetscInt       *nrow,*ncol,*w3,*w4,start;
1563   PetscErrorCode ierr;
1564   PetscMPIInt    size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc;
1565   PetscInt       **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row;
1566   PetscInt       nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol;
1567   PetscInt       **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
1568   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
1569   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
1570   PetscInt       bs     =C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
1571   PetscInt       cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
1572   PetscInt       *bmap  = c->garray,ctmp,rstart=c->rstartbs;
1573   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3;
1574   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3;
1575   MPI_Comm       comm;
1576   PetscBool      flag;
1577   PetscMPIInt    *onodes1,*olengths1;
1578   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
1579 
1580   /* variables below are used for the matrix numerical values - case of !ijonly */
1581   MPI_Request *r_waits4,*s_waits4;
1582   MPI_Status  *r_status4,*s_status4;
1583   MatScalar   **rbuf4,**sbuf_aa,*vals,*mat_a = NULL,*sbuf_aa_i,*vworkA = NULL,*vworkB = NULL;
1584   MatScalar   *a_a=a->a,*b_a=b->a;
1585 
1586 #if defined(PETSC_USE_CTABLE)
1587   PetscInt   tt;
1588   PetscTable *rmap,*cmap,rmap_i,cmap_i=NULL;
1589 #else
1590   PetscInt **cmap,*cmap_i=NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
1591 #endif
1592 
1593   PetscFunctionBegin;
1594   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1595   tag0 = ((PetscObject)C)->tag;
1596   size = c->size;
1597   rank = c->rank;
1598   //if (!rank) printf("MatGetSubMatrices_MPIBAIJ scall %d, bs %d\n",scall,bs);
1599 
1600   /* Get some new tags to keep the communication clean */
1601   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
1602   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
1603   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
1604 
1605 #if defined(PETSC_USE_CTABLE)
1606   ierr = PetscMalloc4(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol);CHKERRQ(ierr);
1607 #else
1608   ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,Mbs+1,&rtable);CHKERRQ(ierr);
1609   /* Create hash table for the mapping :row -> proc*/
1610   for (i=0,j=0; i<size; i++) {
1611     jmax = C->rmap->range[i+1]/bs;
1612     for (; j<jmax; j++) rtable[j] = i;
1613   }
1614 #endif
1615 
1616   for (i=0; i<ismax; i++) {
1617     if (allrows[i]) {
1618       irow[i] = NULL;
1619       nrow[i] = C->rmap->N/bs;
1620     } else {
1621       ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
1622       ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
1623     }
1624 
1625     if (allcolumns[i]) {
1626       icol[i] = NULL;
1627       ncol[i] = C->cmap->N/bs;
1628     } else {
1629       ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
1630       ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
1631     }
1632   }
1633 
1634   /* evaluate communication - mesg to who,length of mesg,and buffer space
1635      required. Based on this, buffers are allocated, and data copied into them*/
1636   ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr);
1637   for (i=0; i<ismax; i++) {
1638     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
1639     jmax   = nrow[i];
1640     irow_i = irow[i];
1641     for (j=0; j<jmax; j++) {
1642       if (allrows[i]) row = j;
1643       else row = irow_i[j];
1644 
1645 #if defined(PETSC_USE_CTABLE)
1646       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1647 #else
1648       proc = rtable[row];
1649 #endif
1650       w4[proc]++;
1651     }
1652     for (j=0; j<size; j++) {
1653       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
1654     }
1655   }
1656 
1657   nrqs     = 0;              /* no of outgoing messages */
1658   msz      = 0;              /* total mesg length for all proc */
1659   w1[rank] = 0;              /* no mesg sent to intself */
1660   w3[rank] = 0;
1661   for (i=0; i<size; i++) {
1662     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
1663   }
1664   ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/
1665   for (i=0,j=0; i<size; i++) {
1666     if (w1[i]) { pa[j] = i; j++; }
1667   }
1668 
1669   /* Each message would have a header = 1 + 2*(no of IS) + data */
1670   for (i=0; i<nrqs; i++) {
1671     j     = pa[i];
1672     w1[j] += w2[j] + 2* w3[j];
1673     msz   += w1[j];
1674   }
1675 
1676   /* Determine the number of messages to expect, their lengths, from from-ids */
1677   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
1678   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
1679 
1680   /* Now post the Irecvs corresponding to these messages */
1681   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
1682 
1683   ierr = PetscFree(onodes1);CHKERRQ(ierr);
1684   ierr = PetscFree(olengths1);CHKERRQ(ierr);
1685 
1686   /* Allocate Memory for outgoing messages */
1687   ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr);
1688   ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
1689   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
1690   {
1691     PetscInt *iptr = tmp,ict = 0;
1692     for (i=0; i<nrqs; i++) {
1693       j        = pa[i];
1694       iptr    += ict;
1695       sbuf1[j] = iptr;
1696       ict      = w1[j];
1697     }
1698   }
1699 
1700   /* Form the outgoing messages */
1701   /* Initialise the header space */
1702   for (i=0; i<nrqs; i++) {
1703     j           = pa[i];
1704     sbuf1[j][0] = 0;
1705     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
1706     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
1707   }
1708 
1709   /* Parse the isrow and copy data into outbuf */
1710   for (i=0; i<ismax; i++) {
1711     ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
1712     irow_i = irow[i];
1713     jmax   = nrow[i];
1714     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
1715       if (allrows[i]) row = j;
1716       else row = irow_i[j];
1717 
1718 #if defined(PETSC_USE_CTABLE)
1719       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1720 #else
1721       proc = rtable[row];
1722 #endif
1723       if (proc != rank) { /* copy to the outgoing buf*/
1724         ctr[proc]++;
1725         *ptr[proc] = row;
1726         ptr[proc]++;
1727       }
1728     }
1729     /* Update the headers for the current IS */
1730     for (j=0; j<size; j++) { /* Can Optimise this loop too */
1731       if ((ctr_j = ctr[j])) {
1732         sbuf1_j        = sbuf1[j];
1733         k              = ++sbuf1_j[0];
1734         sbuf1_j[2*k]   = ctr_j;
1735         sbuf1_j[2*k-1] = i;
1736       }
1737     }
1738   }
1739 
1740   /*  Now  post the sends */
1741   ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr);
1742   for (i=0; i<nrqs; ++i) {
1743     j    = pa[i];
1744     ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
1745   }
1746 
1747   /* Post Recieves to capture the buffer size */
1748   ierr     = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr);
1749   ierr     = PetscMalloc1(nrqs+1,&rbuf2);CHKERRQ(ierr);
1750   rbuf2[0] = tmp + msz;
1751   for (i=1; i<nrqs; ++i) {
1752     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
1753   }
1754   for (i=0; i<nrqs; ++i) {
1755     j        = pa[i];
1756     ierr     = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
1757   }
1758 
1759   /* Send to other procs the buf size they should allocate */
1760 
1761   /* Receive messages*/
1762   ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr);
1763   ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr);
1764   ierr = PetscMalloc3(nrqr+1,&sbuf2,nrqr,&req_size,nrqr,&req_source);CHKERRQ(ierr);
1765   {
1766     Mat_SeqBAIJ *sA  = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
1767     PetscInt    *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
1768 
1769     for (i=0; i<nrqr; ++i) {
1770       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
1771 
1772       req_size[idex] = 0;
1773       rbuf1_i        = rbuf1[idex];
1774       start          = 2*rbuf1_i[0] + 1;
1775       ierr           = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
1776       ierr           = PetscMalloc1(end,&sbuf2[idex]);CHKERRQ(ierr);
1777       sbuf2_i        = sbuf2[idex];
1778       for (j=start; j<end; j++) {
1779         id              = rbuf1_i[j] - rstart;
1780         ncols           = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
1781         sbuf2_i[j]      = ncols;
1782         req_size[idex] += ncols;
1783       }
1784       req_source[idex] = r_status1[i].MPI_SOURCE;
1785       /* form the header */
1786       sbuf2_i[0] = req_size[idex];
1787       for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
1788       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
1789     }
1790   }
1791   ierr = PetscFree(r_status1);CHKERRQ(ierr);
1792   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
1793 
1794   /*  recv buffer sizes */
1795   /* Receive messages*/
1796   ierr = PetscMalloc1(nrqs+1,&rbuf3);CHKERRQ(ierr);
1797   ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr);
1798   ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr);
1799   if (!ijonly) {
1800     ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr);
1801     ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr);
1802   }
1803 
1804   for (i=0; i<nrqs; ++i) {
1805     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
1806     ierr = PetscMalloc1(rbuf2[idex][0],&rbuf3[idex]);CHKERRQ(ierr);
1807     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
1808     if (!ijonly) {
1809       ierr = PetscMalloc1(rbuf2[idex][0]*bs2,&rbuf4[idex]);CHKERRQ(ierr);
1810       ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
1811     }
1812   }
1813   ierr = PetscFree(r_status2);CHKERRQ(ierr);
1814   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
1815 
1816   /* Wait on sends1 and sends2 */
1817   ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr);
1818   ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr);
1819 
1820   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
1821   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
1822   ierr = PetscFree(s_status1);CHKERRQ(ierr);
1823   ierr = PetscFree(s_status2);CHKERRQ(ierr);
1824   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
1825   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
1826 
1827   /* Now allocate buffers for a->j, and send them off */
1828   ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr);
1829   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1830   ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr);
1831   for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
1832 
1833   ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr);
1834   {
1835     for (i=0; i<nrqr; i++) {
1836       rbuf1_i   = rbuf1[i];
1837       sbuf_aj_i = sbuf_aj[i];
1838       ct1       = 2*rbuf1_i[0] + 1;
1839       ct2       = 0;
1840       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
1841         kmax = rbuf1[i][2*j];
1842         for (k=0; k<kmax; k++,ct1++) {
1843           row    = rbuf1_i[ct1] - rstart;
1844           nzA    = a_i[row+1] - a_i[row];
1845           nzB    = b_i[row+1] - b_i[row];
1846           ncols  = nzA + nzB;
1847           cworkA = a_j + a_i[row];
1848           cworkB = b_j + b_i[row];
1849 
1850           /* load the column indices for this row into cols*/
1851           cols = sbuf_aj_i + ct2;
1852           for (l=0; l<nzB; l++) {
1853             if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
1854             else break;
1855           }
1856           imark = l;
1857           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
1858           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
1859           ct2 += ncols;
1860         }
1861       }
1862       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
1863     }
1864   }
1865   ierr = PetscMalloc1(nrqs+1,&r_status3);CHKERRQ(ierr);
1866   ierr = PetscMalloc1(nrqr+1,&s_status3);CHKERRQ(ierr);
1867 
1868   /* Allocate buffers for a->a, and send them off */
1869   if (!ijonly) {
1870     ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr);
1871     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1872     ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr);
1873     for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
1874 
1875     ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr);
1876     {
1877       for (i=0; i<nrqr; i++) {
1878         rbuf1_i   = rbuf1[i];
1879         sbuf_aa_i = sbuf_aa[i];
1880         ct1       = 2*rbuf1_i[0]+1;
1881         ct2       = 0;
1882         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
1883           kmax = rbuf1_i[2*j];
1884           for (k=0; k<kmax; k++,ct1++) {
1885             row    = rbuf1_i[ct1] - rstart;
1886             nzA    = a_i[row+1] - a_i[row];
1887             nzB    = b_i[row+1] - b_i[row];
1888             ncols  = nzA + nzB;
1889             cworkB = b_j + b_i[row];
1890             vworkA = a_a + a_i[row]*bs2;
1891             vworkB = b_a + b_i[row]*bs2;
1892 
1893             /* load the column values for this row into vals*/
1894             vals = sbuf_aa_i+ct2*bs2;
1895             for (l=0; l<nzB; l++) {
1896               if ((bmap[cworkB[l]]) < cstart) {
1897                 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1898               } else break;
1899             }
1900             imark = l;
1901             for (l=0; l<nzA; l++) {
1902               ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1903             }
1904             for (l=imark; l<nzB; l++) {
1905               ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1906             }
1907             ct2 += ncols;
1908           }
1909         }
1910         ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
1911       }
1912     }
1913     ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr);
1914     ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr);
1915   }
1916   ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr);
1917   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
1918 
1919   /* Form the matrix */
1920   /* create col map: global col of C -> local col of submatrices */
1921   {
1922     const PetscInt *icol_i;
1923 #if defined(PETSC_USE_CTABLE)
1924     ierr = PetscMalloc1(1+ismax,&cmap);CHKERRQ(ierr);
1925     for (i=0; i<ismax; i++) {
1926       if (!allcolumns[i]) {
1927         ierr   = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr);
1928         jmax   = ncol[i];
1929         icol_i = icol[i];
1930         cmap_i = cmap[i];
1931         for (j=0; j<jmax; j++) {
1932           ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1933         }
1934       } else {
1935         cmap[i] = NULL;
1936       }
1937     }
1938 #else
1939     ierr = PetscMalloc1(ismax,&cmap);CHKERRQ(ierr);
1940     for (i=0; i<ismax; i++) {
1941       if (!allcolumns[i]) {
1942         ierr   = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr);
1943         jmax   = ncol[i];
1944         icol_i = icol[i];
1945         cmap_i = cmap[i];
1946         for (j=0; j<jmax; j++) {
1947           cmap_i[icol_i[j]] = j+1;
1948         }
1949       } else { /* allcolumns[i] */
1950         cmap[i] = NULL;
1951       }
1952     }
1953 #endif
1954   }
1955 
1956   /* Create lens which is required for MatCreate... */
1957   for (i=0,j=0; i<ismax; i++) j += nrow[i];
1958   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr);
1959   lens[0] = (PetscInt*)(lens + ismax);
1960   ierr    = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr);
1961   for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];
1962 
1963   /* Update lens from local data */
1964   for (i=0; i<ismax; i++) {
1965     jmax = nrow[i];
1966     if (!allcolumns[i]) cmap_i = cmap[i];
1967     irow_i = irow[i];
1968     lens_i = lens[i];
1969     for (j=0; j<jmax; j++) {
1970       if (allrows[i]) row = j;
1971       else row = irow_i[j];
1972 
1973 #if defined(PETSC_USE_CTABLE)
1974       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1975 #else
1976       proc = rtable[row];
1977 #endif
1978       if (proc == rank) {
1979         /* Get indices from matA and then from matB */
1980         row    = row - rstart;
1981         nzA    = a_i[row+1] - a_i[row];
1982         nzB    = b_i[row+1] - b_i[row];
1983         cworkA =  a_j + a_i[row];
1984         cworkB = b_j + b_i[row];
1985         if (!allcolumns[i]) {
1986 #if defined(PETSC_USE_CTABLE)
1987           for (k=0; k<nzA; k++) {
1988             ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1989             if (tt) lens_i[j]++;
1990           }
1991           for (k=0; k<nzB; k++) {
1992             ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1993             if (tt) lens_i[j]++;
1994           }
1995 
1996 #else
1997           for (k=0; k<nzA; k++) {
1998             if (cmap_i[cstart + cworkA[k]]) lens_i[j]++;
1999           }
2000           for (k=0; k<nzB; k++) {
2001             if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++;
2002           }
2003 #endif
2004         } else { /* allcolumns */
2005           lens_i[j] = nzA + nzB;
2006         }
2007       }
2008     }
2009   }
2010 #if defined(PETSC_USE_CTABLE)
2011   /* Create row map*/
2012   ierr = PetscMalloc1(1+ismax,&rmap);CHKERRQ(ierr);
2013   for (i=0; i<ismax; i++) {
2014     ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr);
2015   }
2016 #else
2017   /* Create row map*/
2018   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr);
2019   rmap[0] = (PetscInt*)(rmap + ismax);
2020   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr);
2021   for (i=1; i<ismax; i++) rmap[i] = rmap[i-1] + Mbs;
2022 #endif
2023   for (i=0; i<ismax; i++) {
2024     irow_i = irow[i];
2025     jmax   = nrow[i];
2026 #if defined(PETSC_USE_CTABLE)
2027     rmap_i = rmap[i];
2028     for (j=0; j<jmax; j++) {
2029       if (allrows[i]) {
2030         ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
2031       } else {
2032         ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
2033       }
2034     }
2035 #else
2036     rmap_i = rmap[i];
2037     for (j=0; j<jmax; j++) {
2038       if (allrows[i]) rmap_i[j] = j;
2039       else rmap_i[irow_i[j]] = j;
2040     }
2041 #endif
2042   }
2043 
2044   /* Update lens from offproc data */
2045   {
2046     PetscInt    *rbuf2_i,*rbuf3_i,*sbuf1_i;
2047     PetscMPIInt ii;
2048 
2049     for (tmp2=0; tmp2<nrqs; tmp2++) {
2050       ierr    = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr);
2051       idex    = pa[ii];
2052       sbuf1_i = sbuf1[idex];
2053       jmax    = sbuf1_i[0];
2054       ct1     = 2*jmax+1;
2055       ct2     = 0;
2056       rbuf2_i = rbuf2[ii];
2057       rbuf3_i = rbuf3[ii];
2058       for (j=1; j<=jmax; j++) {
2059         is_no  = sbuf1_i[2*j-1];
2060         max1   = sbuf1_i[2*j];
2061         lens_i = lens[is_no];
2062         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
2063         rmap_i = rmap[is_no];
2064         for (k=0; k<max1; k++,ct1++) {
2065 #if defined(PETSC_USE_CTABLE)
2066           ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
2067           row--;
2068           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
2069 #else
2070           row = rmap_i[sbuf1_i[ct1]];  /* the val in the new matrix to be */
2071 #endif
2072           max2 = rbuf2_i[ct1];
2073           for (l=0; l<max2; l++,ct2++) {
2074             if (!allcolumns[is_no]) {
2075 #if defined(PETSC_USE_CTABLE)
2076               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
2077               if (tt) lens_i[row]++;
2078 #else
2079               if (cmap_i[rbuf3_i[ct2]]) lens_i[row]++;
2080 #endif
2081             } else { /* allcolumns */
2082               lens_i[row]++;
2083             }
2084           }
2085         }
2086       }
2087     }
2088   }
2089   ierr = PetscFree(r_status3);CHKERRQ(ierr);
2090   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
2091   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
2092   ierr = PetscFree(s_status3);CHKERRQ(ierr);
2093   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
2094 
2095   /* Create the submatrices */
2096   if (scall == MAT_REUSE_MATRIX) {
2097     if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet");
2098 
2099     for (i=0; i<ismax; i++) {
2100       mat = (Mat_SeqBAIJ*)(submats[i]->data);
2101       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");
2102       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
2103       if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
2104       /* Initial matrix as if empty */
2105       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr);
2106 
2107       submats[i]->factortype = C->factortype;
2108     }
2109   } else {
2110     PetscInt bs_tmp;
2111     if (ijonly) bs_tmp = 1;
2112     else        bs_tmp = bs;
2113     for (i=0; i<ismax; i++) {
2114       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
2115       ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr);
2116       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
2117       ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr);
2118       ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */
2119     }
2120   }
2121 
2122   /* Assemble the matrices */
2123   /* First assemble the local rows */
2124   {
2125     PetscInt  ilen_row,*imat_ilen,*imat_j,*imat_i;
2126     MatScalar *imat_a = NULL;
2127 
2128     for (i=0; i<ismax; i++) {
2129       mat       = (Mat_SeqBAIJ*)submats[i]->data;
2130       imat_ilen = mat->ilen;
2131       imat_j    = mat->j;
2132       imat_i    = mat->i;
2133       if (!ijonly) imat_a = mat->a;
2134       if (!allcolumns[i]) cmap_i = cmap[i];
2135       rmap_i = rmap[i];
2136       irow_i = irow[i];
2137       jmax   = nrow[i];
2138       for (j=0; j<jmax; j++) {
2139         if (allrows[i]) row = j;
2140         else row = irow_i[j];
2141 
2142 #if defined(PETSC_USE_CTABLE)
2143         ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
2144 #else
2145         proc = rtable[row];
2146 #endif
2147         if (proc == rank) {
2148           row    = row - rstart;
2149           nzA    = a_i[row+1] - a_i[row];
2150           nzB    = b_i[row+1] - b_i[row];
2151           cworkA = a_j + a_i[row];
2152           cworkB = b_j + b_i[row];
2153           if (!ijonly) {
2154             vworkA = a_a + a_i[row]*bs2;
2155             vworkB = b_a + b_i[row]*bs2;
2156           }
2157 #if defined(PETSC_USE_CTABLE)
2158           ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr);
2159           row--;
2160           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
2161 #else
2162           row = rmap_i[row + rstart];
2163 #endif
2164           mat_i = imat_i[row];
2165           if (!ijonly) mat_a = imat_a + mat_i*bs2;
2166           mat_j    = imat_j + mat_i;
2167           ilen_row = imat_ilen[row];
2168 
2169           /* load the column indices for this row into cols*/
2170           if (!allcolumns[i]) {
2171             for (l=0; l<nzB; l++) {
2172               if ((ctmp = bmap[cworkB[l]]) < cstart) {
2173 #if defined(PETSC_USE_CTABLE)
2174                 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr);
2175                 if (tcol) {
2176 #else
2177                 if ((tcol = cmap_i[ctmp])) {
2178 #endif
2179                   *mat_j++ = tcol - 1;
2180                   ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
2181                   mat_a   += bs2;
2182                   ilen_row++;
2183                 }
2184               } else break;
2185             }
2186             imark = l;
2187             for (l=0; l<nzA; l++) {
2188 #if defined(PETSC_USE_CTABLE)
2189               ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
2190               if (tcol) {
2191 #else
2192               if ((tcol = cmap_i[cstart + cworkA[l]])) {
2193 #endif
2194                 *mat_j++ = tcol - 1;
2195                 if (!ijonly) {
2196                   ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
2197                   mat_a += bs2;
2198                 }
2199                 ilen_row++;
2200               }
2201             }
2202             for (l=imark; l<nzB; l++) {
2203 #if defined(PETSC_USE_CTABLE)
2204               ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
2205               if (tcol) {
2206 #else
2207               if ((tcol = cmap_i[bmap[cworkB[l]]])) {
2208 #endif
2209                 *mat_j++ = tcol - 1;
2210                 if (!ijonly) {
2211                   ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
2212                   mat_a += bs2;
2213                 }
2214                 ilen_row++;
2215               }
2216             }
2217           } else { /* allcolumns */
2218             for (l=0; l<nzB; l++) {
2219               if ((ctmp = bmap[cworkB[l]]) < cstart) {
2220                 *mat_j++ = ctmp;
2221                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
2222                 mat_a   += bs2;
2223                 ilen_row++;
2224               } else break;
2225             }
2226             imark = l;
2227             for (l=0; l<nzA; l++) {
2228               *mat_j++ = cstart+cworkA[l];
2229               if (!ijonly) {
2230                 ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
2231                 mat_a += bs2;
2232               }
2233               ilen_row++;
2234             }
2235             for (l=imark; l<nzB; l++) {
2236               *mat_j++ = bmap[cworkB[l]];
2237               if (!ijonly) {
2238                 ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
2239                 mat_a += bs2;
2240               }
2241               ilen_row++;
2242             }
2243           }
2244           imat_ilen[row] = ilen_row;
2245         }
2246       }
2247     }
2248   }
2249 
2250   /*   Now assemble the off proc rows*/
2251   {
2252     PetscInt    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
2253     PetscInt    *imat_j,*imat_i;
2254     MatScalar   *imat_a = NULL,*rbuf4_i = NULL;
2255     PetscMPIInt ii;
2256 
2257     for (tmp2=0; tmp2<nrqs; tmp2++) {
2258       if (ijonly) ii = tmp2;
2259       else {
2260         ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr);
2261       }
2262       idex    = pa[ii];
2263       sbuf1_i = sbuf1[idex];
2264       jmax    = sbuf1_i[0];
2265       ct1     = 2*jmax + 1;
2266       ct2     = 0;
2267       rbuf2_i = rbuf2[ii];
2268       rbuf3_i = rbuf3[ii];
2269       if (!ijonly) rbuf4_i = rbuf4[ii];
2270       for (j=1; j<=jmax; j++) {
2271         is_no = sbuf1_i[2*j-1];
2272         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
2273         rmap_i    = rmap[is_no];
2274         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
2275         imat_ilen = mat->ilen;
2276         imat_j    = mat->j;
2277         imat_i    = mat->i;
2278         if (!ijonly) imat_a = mat->a;
2279         max1 = sbuf1_i[2*j];
2280         for (k=0; k<max1; k++,ct1++) {
2281           row = sbuf1_i[ct1];
2282 #if defined(PETSC_USE_CTABLE)
2283           ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr);
2284           row--;
2285           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
2286 #else
2287           row = rmap_i[row];
2288 #endif
2289           ilen  = imat_ilen[row];
2290           mat_i = imat_i[row];
2291           if (!ijonly) mat_a = imat_a + mat_i*bs2;
2292           mat_j = imat_j + mat_i;
2293           max2  = rbuf2_i[ct1];
2294 
2295           if (!allcolumns[is_no]) {
2296             for (l=0; l<max2; l++,ct2++) {
2297 #if defined(PETSC_USE_CTABLE)
2298               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
2299               if (tcol) {
2300 #else
2301               if ((tcol = cmap_i[rbuf3_i[ct2]])) {
2302 #endif
2303                 *mat_j++ = tcol - 1;
2304                 if (!ijonly) {
2305                   ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
2306                   mat_a += bs2;
2307                 }
2308                 ilen++;
2309               }
2310             }
2311           } else { /* allcolumns */
2312             for (l=0; l<max2; l++,ct2++) {
2313               *mat_j++ = rbuf3_i[ct2];
2314               if (!ijonly) {
2315                 ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
2316                 mat_a += bs2;
2317               }
2318               ilen++;
2319             }
2320           }
2321           imat_ilen[row] = ilen;
2322         }
2323       }
2324     }
2325   }
2326   /* sort the rows */
2327   {
2328     PetscInt  ilen_row,*imat_ilen,*imat_j,*imat_i;
2329     MatScalar *imat_a = NULL;
2330     MatScalar *work;
2331 
2332     ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr);
2333     for (i=0; i<ismax; i++) {
2334       mat       = (Mat_SeqBAIJ*)submats[i]->data;
2335       imat_ilen = mat->ilen;
2336       imat_j    = mat->j;
2337       imat_i    = mat->i;
2338       if (!ijonly) imat_a = mat->a;
2339       if (allcolumns[i]) continue;
2340       jmax   = nrow[i];
2341       for (j=0; j<jmax; j++) {
2342         mat_i = imat_i[j];
2343         if (!ijonly) mat_a = imat_a + mat_i*bs2;
2344         mat_j    = imat_j + mat_i;
2345         ilen_row = imat_ilen[j];
2346         if (!ijonly) {ierr = PetscSortIntWithDataArray(ilen_row,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);}
2347         else {ierr = PetscSortInt(ilen_row,mat_j);CHKERRQ(ierr);}
2348       }
2349     }
2350     ierr = PetscFree(work);CHKERRQ(ierr);
2351   }
2352   if (!ijonly) {
2353     ierr = PetscFree(r_status4);CHKERRQ(ierr);
2354     ierr = PetscFree(r_waits4);CHKERRQ(ierr);
2355     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
2356     ierr = PetscFree(s_waits4);CHKERRQ(ierr);
2357     ierr = PetscFree(s_status4);CHKERRQ(ierr);
2358   }
2359 
2360   /* Restore the indices */
2361   for (i=0; i<ismax; i++) {
2362     if (!allrows[i]) {
2363       ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
2364     }
2365     if (!allcolumns[i]) {
2366       ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
2367     }
2368   }
2369 
2370   /* Destroy allocated memory */
2371 #if defined(PETSC_USE_CTABLE)
2372   ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr);
2373 #else
2374   ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr);
2375 #endif
2376   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
2377   ierr = PetscFree(pa);CHKERRQ(ierr);
2378 
2379   ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
2380   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
2381   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
2382   for (i=0; i<nrqr; ++i) {
2383     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
2384   }
2385   for (i=0; i<nrqs; ++i) {
2386     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
2387   }
2388   ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr);
2389   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
2390   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
2391   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
2392   if (!ijonly) {
2393     for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);}
2394     ierr = PetscFree(rbuf4);CHKERRQ(ierr);
2395     ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
2396     ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
2397   }
2398 
2399 #if defined(PETSC_USE_CTABLE)
2400   for (i=0; i<ismax; i++) {
2401     ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr);
2402   }
2403 #endif
2404   ierr = PetscFree(rmap);CHKERRQ(ierr);
2405 
2406   for (i=0; i<ismax; i++) {
2407     if (!allcolumns[i]) {
2408 #if defined(PETSC_USE_CTABLE)
2409       ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr);
2410 #else
2411       ierr = PetscFree(cmap[i]);CHKERRQ(ierr);
2412 #endif
2413     }
2414   }
2415   ierr = PetscFree(cmap);CHKERRQ(ierr);
2416   ierr = PetscFree(lens);CHKERRQ(ierr);
2417 
2418   for (i=0; i<ismax; i++) {
2419     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2420     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2421   }
2422 
2423   c->ijonly = PETSC_FALSE; /* set back to the default */
2424   PetscFunctionReturn(0);
2425 }
2426 
2427