xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 2f6eced2a19e978d64f27de66fbc6c26cc5c7934)
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_new,*iscol_new;
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_new,ismax+1,&iscol_new);CHKERRQ(ierr);
530   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_new);CHKERRQ(ierr);
531   ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_new);CHKERRQ(ierr);
532 
533   /* Check for special case: each processor gets entire matrix columns */
534   ierr = PetscMalloc2(ismax+1,&allcolumns,ismax+1,&allrows);CHKERRQ(ierr);
535   for (i=0; i<ismax; i++) {
536     ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr);
537     ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr);
538     if (colflag && ncol == C->cmap->N) allcolumns[i] = PETSC_TRUE;
539     else allcolumns[i] = PETSC_FALSE;
540 
541     ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr);
542     ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr);
543     if (colflag && nrow == C->rmap->N) allrows[i] = PETSC_TRUE;
544     else allrows[i] = PETSC_FALSE;
545   }
546 
547   /* Allocate memory to hold all the submatrices */
548   if (scall != MAT_REUSE_MATRIX) {
549     ierr = PetscMalloc1(ismax+1,submat);CHKERRQ(ierr);
550   }
551   /* Determine the number of stages through which submatrices are done */
552   nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt));
553   if (!nmax) nmax = 1;
554   nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0);
555 
556   /* Make sure every processor loops through the nstages */
557   ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr);
558   for (i=0,pos=0; i<nstages; i++) {
559     if (pos+nmax <= ismax) max_no = nmax;
560     else if (pos == ismax) max_no = 0;
561     else                   max_no = ismax-pos;
562     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr);
563     pos += max_no;
564   }
565 
566   for (i=0; i<ismax; i++) {
567     ierr = ISDestroy(&isrow_new[i]);CHKERRQ(ierr);
568     ierr = ISDestroy(&iscol_new[i]);CHKERRQ(ierr);
569   }
570   ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr);
571   ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr);
572   PetscFunctionReturn(0);
573 }
574 
575 #if defined(PETSC_USE_CTABLE)
576 PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
577 {
578   PetscInt       nGlobalNd = proc_gnode[size];
579   PetscMPIInt    fproc;
580   PetscErrorCode ierr;
581 
582   PetscFunctionBegin;
583   ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr);
584   if (fproc > size) fproc = size;
585   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
586     if (row < proc_gnode[fproc]) fproc--;
587     else                         fproc++;
588   }
589   *rank = fproc;
590   PetscFunctionReturn(0);
591 }
592 #endif
593 
594 /* -------------------------------------------------------------------------*/
595 /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
596 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats)
597 {
598   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
599   Mat            A  = c->A;
600   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
601   const PetscInt **irow,**icol,*irow_i;
602   PetscInt       *nrow,*ncol,*w3,*w4,start;
603   PetscErrorCode ierr;
604   PetscMPIInt    size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc;
605   PetscInt       **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row;
606   PetscInt       nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol;
607   PetscInt       **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
608   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
609   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
610   PetscInt       bs     =C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
611   PetscInt       cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
612   PetscInt       *bmap  = c->garray,ctmp,rstart=c->rstartbs;
613   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3;
614   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3;
615   MPI_Comm       comm;
616   PetscBool      flag;
617   PetscMPIInt    *onodes1,*olengths1;
618   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
619 
620   /* variables below are used for the matrix numerical values - case of !ijonly */
621   MPI_Request *r_waits4,*s_waits4;
622   MPI_Status  *r_status4,*s_status4;
623   MatScalar   **rbuf4,**sbuf_aa,*vals,*mat_a = NULL,*sbuf_aa_i,*vworkA = NULL,*vworkB = NULL;
624   MatScalar   *a_a=a->a,*b_a=b->a;
625 
626 #if defined(PETSC_USE_CTABLE)
627   PetscInt   tt;
628   PetscTable *rmap,*cmap,rmap_i,cmap_i=NULL;
629 #else
630   PetscInt **cmap,*cmap_i=NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
631 #endif
632 
633   PetscFunctionBegin;
634   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
635   tag0 = ((PetscObject)C)->tag;
636   size = c->size;
637   rank = c->rank;
638 
639   /* Get some new tags to keep the communication clean */
640   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
641   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
642   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
643 
644 #if defined(PETSC_USE_CTABLE)
645   ierr = PetscMalloc4(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol);CHKERRQ(ierr);
646 #else
647   ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,Mbs+1,&rtable);CHKERRQ(ierr);
648   /* Create hash table for the mapping :row -> proc*/
649   for (i=0,j=0; i<size; i++) {
650     jmax = C->rmap->range[i+1]/bs;
651     for (; j<jmax; j++) rtable[j] = i;
652   }
653 #endif
654 
655   for (i=0; i<ismax; i++) {
656     if (allrows[i]) {
657       irow[i] = NULL;
658       nrow[i] = C->rmap->N/bs;
659     } else {
660       ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
661       ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
662     }
663 
664     if (allcolumns[i]) {
665       icol[i] = NULL;
666       ncol[i] = C->cmap->N/bs;
667     } else {
668       ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
669       ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
670     }
671   }
672 
673   /* evaluate communication - mesg to who,length of mesg,and buffer space
674      required. Based on this, buffers are allocated, and data copied into them*/
675   ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr);
676   for (i=0; i<ismax; i++) {
677     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
678     jmax   = nrow[i];
679     irow_i = irow[i];
680     for (j=0; j<jmax; j++) {
681       if (allrows[i]) row = j;
682       else row = irow_i[j];
683 
684 #if defined(PETSC_USE_CTABLE)
685       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
686 #else
687       proc = rtable[row];
688 #endif
689       w4[proc]++;
690     }
691     for (j=0; j<size; j++) {
692       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
693     }
694   }
695 
696   nrqs     = 0;              /* no of outgoing messages */
697   msz      = 0;              /* total mesg length for all proc */
698   w1[rank] = 0;              /* no mesg sent to intself */
699   w3[rank] = 0;
700   for (i=0; i<size; i++) {
701     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
702   }
703   ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/
704   for (i=0,j=0; i<size; i++) {
705     if (w1[i]) { pa[j] = i; j++; }
706   }
707 
708   /* Each message would have a header = 1 + 2*(no of IS) + data */
709   for (i=0; i<nrqs; i++) {
710     j     = pa[i];
711     w1[j] += w2[j] + 2* w3[j];
712     msz   += w1[j];
713   }
714 
715   /* Determine the number of messages to expect, their lengths, from from-ids */
716   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
717   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
718 
719   /* Now post the Irecvs corresponding to these messages */
720   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
721 
722   ierr = PetscFree(onodes1);CHKERRQ(ierr);
723   ierr = PetscFree(olengths1);CHKERRQ(ierr);
724 
725   /* Allocate Memory for outgoing messages */
726   ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr);
727   ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
728   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
729   {
730     PetscInt *iptr = tmp,ict = 0;
731     for (i=0; i<nrqs; i++) {
732       j        = pa[i];
733       iptr    += ict;
734       sbuf1[j] = iptr;
735       ict      = w1[j];
736     }
737   }
738 
739   /* Form the outgoing messages */
740   /* Initialise the header space */
741   for (i=0; i<nrqs; i++) {
742     j           = pa[i];
743     sbuf1[j][0] = 0;
744     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
745     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
746   }
747 
748   /* Parse the isrow and copy data into outbuf */
749   for (i=0; i<ismax; i++) {
750     ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
751     irow_i = irow[i];
752     jmax   = nrow[i];
753     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
754       if (allrows[i]) row = j;
755       else row = irow_i[j];
756 
757 #if defined(PETSC_USE_CTABLE)
758       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
759 #else
760       proc = rtable[row];
761 #endif
762       if (proc != rank) { /* copy to the outgoing buf*/
763         ctr[proc]++;
764         *ptr[proc] = row;
765         ptr[proc]++;
766       }
767     }
768     /* Update the headers for the current IS */
769     for (j=0; j<size; j++) { /* Can Optimise this loop too */
770       if ((ctr_j = ctr[j])) {
771         sbuf1_j        = sbuf1[j];
772         k              = ++sbuf1_j[0];
773         sbuf1_j[2*k]   = ctr_j;
774         sbuf1_j[2*k-1] = i;
775       }
776     }
777   }
778 
779   /*  Now  post the sends */
780   ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr);
781   for (i=0; i<nrqs; ++i) {
782     j    = pa[i];
783     ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
784   }
785 
786   /* Post Recieves to capture the buffer size */
787   ierr     = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr);
788   ierr     = PetscMalloc1(nrqs+1,&rbuf2);CHKERRQ(ierr);
789   rbuf2[0] = tmp + msz;
790   for (i=1; i<nrqs; ++i) {
791     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
792   }
793   for (i=0; i<nrqs; ++i) {
794     j        = pa[i];
795     ierr     = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
796   }
797 
798   /* Send to other procs the buf size they should allocate */
799 
800   /* Receive messages*/
801   ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr);
802   ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr);
803   ierr = PetscMalloc3(nrqr+1,&sbuf2,nrqr,&req_size,nrqr,&req_source);CHKERRQ(ierr);
804   {
805     Mat_SeqBAIJ *sA  = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
806     PetscInt    *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
807 
808     for (i=0; i<nrqr; ++i) {
809       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
810 
811       req_size[idex] = 0;
812       rbuf1_i        = rbuf1[idex];
813       start          = 2*rbuf1_i[0] + 1;
814       ierr           = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
815       ierr           = PetscMalloc1(end,&sbuf2[idex]);CHKERRQ(ierr);
816       sbuf2_i        = sbuf2[idex];
817       for (j=start; j<end; j++) {
818         id              = rbuf1_i[j] - rstart;
819         ncols           = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
820         sbuf2_i[j]      = ncols;
821         req_size[idex] += ncols;
822       }
823       req_source[idex] = r_status1[i].MPI_SOURCE;
824       /* form the header */
825       sbuf2_i[0] = req_size[idex];
826       for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
827       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
828     }
829   }
830   ierr = PetscFree(r_status1);CHKERRQ(ierr);
831   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
832 
833   /*  recv buffer sizes */
834   /* Receive messages*/
835   ierr = PetscMalloc1(nrqs+1,&rbuf3);CHKERRQ(ierr);
836   ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr);
837   ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr);
838   if (!ijonly) {
839     ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr);
840     ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr);
841   }
842 
843   for (i=0; i<nrqs; ++i) {
844     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
845     ierr = PetscMalloc1(rbuf2[idex][0],&rbuf3[idex]);CHKERRQ(ierr);
846     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
847     if (!ijonly) {
848       ierr = PetscMalloc1(rbuf2[idex][0]*bs2,&rbuf4[idex]);CHKERRQ(ierr);
849       ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
850     }
851   }
852   ierr = PetscFree(r_status2);CHKERRQ(ierr);
853   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
854 
855   /* Wait on sends1 and sends2 */
856   ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr);
857   ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr);
858 
859   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
860   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
861   ierr = PetscFree(s_status1);CHKERRQ(ierr);
862   ierr = PetscFree(s_status2);CHKERRQ(ierr);
863   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
864   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
865 
866   /* Now allocate buffers for a->j, and send them off */
867   ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr);
868   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
869   ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr);
870   for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
871 
872   ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr);
873   {
874     for (i=0; i<nrqr; i++) {
875       rbuf1_i   = rbuf1[i];
876       sbuf_aj_i = sbuf_aj[i];
877       ct1       = 2*rbuf1_i[0] + 1;
878       ct2       = 0;
879       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
880         kmax = rbuf1[i][2*j];
881         for (k=0; k<kmax; k++,ct1++) {
882           row    = rbuf1_i[ct1] - rstart;
883           nzA    = a_i[row+1] - a_i[row];
884           nzB    = b_i[row+1] - b_i[row];
885           ncols  = nzA + nzB;
886           cworkA = a_j + a_i[row];
887           cworkB = b_j + b_i[row];
888 
889           /* load the column indices for this row into cols*/
890           cols = sbuf_aj_i + ct2;
891           for (l=0; l<nzB; l++) {
892             if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
893             else break;
894           }
895           imark = l;
896           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
897           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
898           ct2 += ncols;
899         }
900       }
901       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
902     }
903   }
904   ierr = PetscMalloc1(nrqs+1,&r_status3);CHKERRQ(ierr);
905   ierr = PetscMalloc1(nrqr+1,&s_status3);CHKERRQ(ierr);
906 
907   /* Allocate buffers for a->a, and send them off */
908   if (!ijonly) {
909     ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr);
910     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
911     ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr);
912     for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
913 
914     ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr);
915     {
916       for (i=0; i<nrqr; i++) {
917         rbuf1_i   = rbuf1[i];
918         sbuf_aa_i = sbuf_aa[i];
919         ct1       = 2*rbuf1_i[0]+1;
920         ct2       = 0;
921         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
922           kmax = rbuf1_i[2*j];
923           for (k=0; k<kmax; k++,ct1++) {
924             row    = rbuf1_i[ct1] - rstart;
925             nzA    = a_i[row+1] - a_i[row];
926             nzB    = b_i[row+1] - b_i[row];
927             ncols  = nzA + nzB;
928             cworkB = b_j + b_i[row];
929             vworkA = a_a + a_i[row]*bs2;
930             vworkB = b_a + b_i[row]*bs2;
931 
932             /* load the column values for this row into vals*/
933             vals = sbuf_aa_i+ct2*bs2;
934             for (l=0; l<nzB; l++) {
935               if ((bmap[cworkB[l]]) < cstart) {
936                 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
937               } else break;
938             }
939             imark = l;
940             for (l=0; l<nzA; l++) {
941               ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
942             }
943             for (l=imark; l<nzB; l++) {
944               ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
945             }
946             ct2 += ncols;
947           }
948         }
949         ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
950       }
951     }
952     ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr);
953     ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr);
954   }
955   ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr);
956   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
957 
958   /* Form the matrix */
959   /* create col map: global col of C -> local col of submatrices */
960   {
961     const PetscInt *icol_i;
962 #if defined(PETSC_USE_CTABLE)
963     ierr = PetscMalloc1(1+ismax,&cmap);CHKERRQ(ierr);
964     for (i=0; i<ismax; i++) {
965       if (!allcolumns[i]) {
966         ierr   = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr);
967         jmax   = ncol[i];
968         icol_i = icol[i];
969         cmap_i = cmap[i];
970         for (j=0; j<jmax; j++) {
971           ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
972         }
973       } else {
974         cmap[i] = NULL;
975       }
976     }
977 #else
978     ierr = PetscMalloc1(ismax,&cmap);CHKERRQ(ierr);
979     for (i=0; i<ismax; i++) {
980       if (!allcolumns[i]) {
981         ierr   = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr);
982         jmax   = ncol[i];
983         icol_i = icol[i];
984         cmap_i = cmap[i];
985         for (j=0; j<jmax; j++) {
986           cmap_i[icol_i[j]] = j+1;
987         }
988       } else { /* allcolumns[i] */
989         cmap[i] = NULL;
990       }
991     }
992 #endif
993   }
994 
995   /* Create lens which is required for MatCreate... */
996   for (i=0,j=0; i<ismax; i++) j += nrow[i];
997   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr);
998   lens[0] = (PetscInt*)(lens + ismax);
999   ierr    = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr);
1000   for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];
1001 
1002   /* Update lens from local data */
1003   for (i=0; i<ismax; i++) {
1004     jmax = nrow[i];
1005     if (!allcolumns[i]) cmap_i = cmap[i];
1006     irow_i = irow[i];
1007     lens_i = lens[i];
1008     for (j=0; j<jmax; j++) {
1009       if (allrows[i]) row = j;
1010       else row = irow_i[j];
1011 
1012 #if defined(PETSC_USE_CTABLE)
1013       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1014 #else
1015       proc = rtable[row];
1016 #endif
1017       if (proc == rank) {
1018         /* Get indices from matA and then from matB */
1019         row    = row - rstart;
1020         nzA    = a_i[row+1] - a_i[row];
1021         nzB    = b_i[row+1] - b_i[row];
1022         cworkA =  a_j + a_i[row];
1023         cworkB = b_j + b_i[row];
1024         if (!allcolumns[i]) {
1025 #if defined(PETSC_USE_CTABLE)
1026           for (k=0; k<nzA; k++) {
1027             ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1028             if (tt) lens_i[j]++;
1029           }
1030           for (k=0; k<nzB; k++) {
1031             ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1032             if (tt) lens_i[j]++;
1033           }
1034 
1035 #else
1036           for (k=0; k<nzA; k++) {
1037             if (cmap_i[cstart + cworkA[k]]) lens_i[j]++;
1038           }
1039           for (k=0; k<nzB; k++) {
1040             if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++;
1041           }
1042 #endif
1043         } else { /* allcolumns */
1044           lens_i[j] = nzA + nzB;
1045         }
1046       }
1047     }
1048   }
1049 #if defined(PETSC_USE_CTABLE)
1050   /* Create row map*/
1051   ierr = PetscMalloc1(1+ismax,&rmap);CHKERRQ(ierr);
1052   for (i=0; i<ismax; i++) {
1053     ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr);
1054   }
1055 #else
1056   /* Create row map*/
1057   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr);
1058   rmap[0] = (PetscInt*)(rmap + ismax);
1059   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr);
1060   for (i=1; i<ismax; i++) rmap[i] = rmap[i-1] + Mbs;
1061 #endif
1062   for (i=0; i<ismax; i++) {
1063     irow_i = irow[i];
1064     jmax   = nrow[i];
1065 #if defined(PETSC_USE_CTABLE)
1066     rmap_i = rmap[i];
1067     for (j=0; j<jmax; j++) {
1068       if (allrows[i]) {
1069         ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1070       } else {
1071         ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1072       }
1073     }
1074 #else
1075     rmap_i = rmap[i];
1076     for (j=0; j<jmax; j++) {
1077       if (allrows[i]) rmap_i[j] = j;
1078       else rmap_i[irow_i[j]] = j;
1079     }
1080 #endif
1081   }
1082 
1083   /* Update lens from offproc data */
1084   {
1085     PetscInt    *rbuf2_i,*rbuf3_i,*sbuf1_i;
1086     PetscMPIInt ii;
1087 
1088     for (tmp2=0; tmp2<nrqs; tmp2++) {
1089       ierr    = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr);
1090       idex    = pa[ii];
1091       sbuf1_i = sbuf1[idex];
1092       jmax    = sbuf1_i[0];
1093       ct1     = 2*jmax+1;
1094       ct2     = 0;
1095       rbuf2_i = rbuf2[ii];
1096       rbuf3_i = rbuf3[ii];
1097       for (j=1; j<=jmax; j++) {
1098         is_no  = sbuf1_i[2*j-1];
1099         max1   = sbuf1_i[2*j];
1100         lens_i = lens[is_no];
1101         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1102         rmap_i = rmap[is_no];
1103         for (k=0; k<max1; k++,ct1++) {
1104 #if defined(PETSC_USE_CTABLE)
1105           ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1106           row--;
1107           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1108 #else
1109           row = rmap_i[sbuf1_i[ct1]];  /* the val in the new matrix to be */
1110 #endif
1111           max2 = rbuf2_i[ct1];
1112           for (l=0; l<max2; l++,ct2++) {
1113             if (!allcolumns[is_no]) {
1114 #if defined(PETSC_USE_CTABLE)
1115               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
1116               if (tt) lens_i[row]++;
1117 #else
1118               if (cmap_i[rbuf3_i[ct2]]) lens_i[row]++;
1119 #endif
1120             } else { /* allcolumns */
1121               lens_i[row]++;
1122             }
1123           }
1124         }
1125       }
1126     }
1127   }
1128   ierr = PetscFree(r_status3);CHKERRQ(ierr);
1129   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
1130   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
1131   ierr = PetscFree(s_status3);CHKERRQ(ierr);
1132   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1133 
1134   /* Create the submatrices */
1135   if (scall == MAT_REUSE_MATRIX) {
1136     if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet");
1137     /*
1138         Assumes new rows are same length as the old rows, hence bug!
1139     */
1140     for (i=0; i<ismax; i++) {
1141       mat = (Mat_SeqBAIJ*)(submats[i]->data);
1142       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");
1143       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
1144       if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1145       /* Initial matrix as if empty */
1146       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr);
1147 
1148       submats[i]->factortype = C->factortype;
1149     }
1150   } else {
1151     PetscInt bs_tmp;
1152     if (ijonly) bs_tmp = 1;
1153     else        bs_tmp = bs;
1154     for (i=0; i<ismax; i++) {
1155       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
1156       ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr);
1157       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
1158       ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr);
1159       ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */
1160     }
1161   }
1162 
1163   /* Assemble the matrices */
1164   /* First assemble the local rows */
1165   {
1166     PetscInt  ilen_row,*imat_ilen,*imat_j,*imat_i;
1167     MatScalar *imat_a = NULL;
1168 
1169     for (i=0; i<ismax; i++) {
1170       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1171       imat_ilen = mat->ilen;
1172       imat_j    = mat->j;
1173       imat_i    = mat->i;
1174       if (!ijonly) imat_a = mat->a;
1175       if (!allcolumns[i]) cmap_i = cmap[i];
1176       rmap_i = rmap[i];
1177       irow_i = irow[i];
1178       jmax   = nrow[i];
1179       for (j=0; j<jmax; j++) {
1180         if (allrows[i]) row = j;
1181         else row = irow_i[j];
1182 
1183 #if defined(PETSC_USE_CTABLE)
1184         ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1185 #else
1186         proc = rtable[row];
1187 #endif
1188         if (proc == rank) {
1189           row    = row - rstart;
1190           nzA    = a_i[row+1] - a_i[row];
1191           nzB    = b_i[row+1] - b_i[row];
1192           cworkA = a_j + a_i[row];
1193           cworkB = b_j + b_i[row];
1194           if (!ijonly) {
1195             vworkA = a_a + a_i[row]*bs2;
1196             vworkB = b_a + b_i[row]*bs2;
1197           }
1198 #if defined(PETSC_USE_CTABLE)
1199           ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr);
1200           row--;
1201           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1202 #else
1203           row = rmap_i[row + rstart];
1204 #endif
1205           mat_i = imat_i[row];
1206           if (!ijonly) mat_a = imat_a + mat_i*bs2;
1207           mat_j    = imat_j + mat_i;
1208           ilen_row = imat_ilen[row];
1209 
1210           /* load the column indices for this row into cols*/
1211           if (!allcolumns[i]) {
1212             for (l=0; l<nzB; l++) {
1213               if ((ctmp = bmap[cworkB[l]]) < cstart) {
1214 #if defined(PETSC_USE_CTABLE)
1215                 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr);
1216                 if (tcol) {
1217 #else
1218                 if ((tcol = cmap_i[ctmp])) {
1219 #endif
1220                   *mat_j++ = tcol - 1;
1221                   ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1222                   mat_a   += bs2;
1223                   ilen_row++;
1224                 }
1225               } else break;
1226             }
1227             imark = l;
1228             for (l=0; l<nzA; l++) {
1229 #if defined(PETSC_USE_CTABLE)
1230               ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1231               if (tcol) {
1232 #else
1233               if ((tcol = cmap_i[cstart + cworkA[l]])) {
1234 #endif
1235                 *mat_j++ = tcol - 1;
1236                 if (!ijonly) {
1237                   ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1238                   mat_a += bs2;
1239                 }
1240                 ilen_row++;
1241               }
1242             }
1243             for (l=imark; l<nzB; l++) {
1244 #if defined(PETSC_USE_CTABLE)
1245               ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1246               if (tcol) {
1247 #else
1248               if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1249 #endif
1250                 *mat_j++ = tcol - 1;
1251                 if (!ijonly) {
1252                   ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1253                   mat_a += bs2;
1254                 }
1255                 ilen_row++;
1256               }
1257             }
1258           } else { /* allcolumns */
1259             for (l=0; l<nzB; l++) {
1260               if ((ctmp = bmap[cworkB[l]]) < cstart) {
1261                 *mat_j++ = ctmp;
1262                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1263                 mat_a   += bs2;
1264                 ilen_row++;
1265               } else break;
1266             }
1267             imark = l;
1268             for (l=0; l<nzA; l++) {
1269               *mat_j++ = cstart+cworkA[l];
1270               if (!ijonly) {
1271                 ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1272                 mat_a += bs2;
1273               }
1274               ilen_row++;
1275             }
1276             for (l=imark; l<nzB; l++) {
1277               *mat_j++ = bmap[cworkB[l]];
1278               if (!ijonly) {
1279                 ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1280                 mat_a += bs2;
1281               }
1282               ilen_row++;
1283             }
1284           }
1285           imat_ilen[row] = ilen_row;
1286         }
1287       }
1288     }
1289   }
1290 
1291   /*   Now assemble the off proc rows*/
1292   {
1293     PetscInt    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1294     PetscInt    *imat_j,*imat_i;
1295     MatScalar   *imat_a = NULL,*rbuf4_i = NULL;
1296     PetscMPIInt ii;
1297 
1298     for (tmp2=0; tmp2<nrqs; tmp2++) {
1299       if (ijonly) ii = tmp2;
1300       else {
1301         ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr);
1302       }
1303       idex    = pa[ii];
1304       sbuf1_i = sbuf1[idex];
1305       jmax    = sbuf1_i[0];
1306       ct1     = 2*jmax + 1;
1307       ct2     = 0;
1308       rbuf2_i = rbuf2[ii];
1309       rbuf3_i = rbuf3[ii];
1310       if (!ijonly) rbuf4_i = rbuf4[ii];
1311       for (j=1; j<=jmax; j++) {
1312         is_no = sbuf1_i[2*j-1];
1313         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1314         rmap_i    = rmap[is_no];
1315         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1316         imat_ilen = mat->ilen;
1317         imat_j    = mat->j;
1318         imat_i    = mat->i;
1319         if (!ijonly) imat_a = mat->a;
1320         max1 = sbuf1_i[2*j];
1321         for (k=0; k<max1; k++,ct1++) {
1322           row = sbuf1_i[ct1];
1323 #if defined(PETSC_USE_CTABLE)
1324           ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr);
1325           row--;
1326           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1327 #else
1328           row = rmap_i[row];
1329 #endif
1330           ilen  = imat_ilen[row];
1331           mat_i = imat_i[row];
1332           if (!ijonly) mat_a = imat_a + mat_i*bs2;
1333           mat_j = imat_j + mat_i;
1334           max2  = rbuf2_i[ct1];
1335 
1336           if (!allcolumns[is_no]) {
1337             for (l=0; l<max2; l++,ct2++) {
1338 #if defined(PETSC_USE_CTABLE)
1339               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1340               if (tcol) {
1341 #else
1342               if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1343 #endif
1344                 *mat_j++ = tcol - 1;
1345                 if (!ijonly) {
1346                   ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1347                   mat_a += bs2;
1348                 }
1349                 ilen++;
1350               }
1351             }
1352           } else { /* allcolumns */
1353             for (l=0; l<max2; l++,ct2++) {
1354               *mat_j++ = rbuf3_i[ct2];
1355               if (!ijonly) {
1356                 ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1357                 mat_a += bs2;
1358               }
1359               ilen++;
1360             }
1361           }
1362           imat_ilen[row] = ilen;
1363         }
1364       }
1365     }
1366   }
1367   /* sort the rows */
1368   {
1369     PetscInt  ilen_row,*imat_ilen,*imat_j,*imat_i;
1370     MatScalar *imat_a = NULL;
1371     MatScalar *work;
1372 
1373     ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr);
1374     for (i=0; i<ismax; i++) {
1375       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1376       imat_ilen = mat->ilen;
1377       imat_j    = mat->j;
1378       imat_i    = mat->i;
1379       if (!ijonly) imat_a = mat->a;
1380       if (allcolumns[i]) continue;
1381       jmax   = nrow[i];
1382       for (j=0; j<jmax; j++) {
1383         mat_i = imat_i[j];
1384         if (!ijonly) mat_a = imat_a + mat_i*bs2;
1385         mat_j    = imat_j + mat_i;
1386         ilen_row = imat_ilen[j];
1387         if (!ijonly) {ierr = PetscSortIntWithDataArray(ilen_row,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);}
1388         else {ierr = PetscSortInt(ilen_row,mat_j);CHKERRQ(ierr);}
1389       }
1390     }
1391     ierr = PetscFree(work);CHKERRQ(ierr);
1392   }
1393   if (!ijonly) {
1394     ierr = PetscFree(r_status4);CHKERRQ(ierr);
1395     ierr = PetscFree(r_waits4);CHKERRQ(ierr);
1396     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
1397     ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1398     ierr = PetscFree(s_status4);CHKERRQ(ierr);
1399   }
1400 
1401   /* Restore the indices */
1402   for (i=0; i<ismax; i++) {
1403     if (!allrows[i]) {
1404       ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1405     }
1406     if (!allcolumns[i]) {
1407       ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1408     }
1409   }
1410 
1411   /* Destroy allocated memory */
1412 #if defined(PETSC_USE_CTABLE)
1413   ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr);
1414 #else
1415   ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr);
1416 #endif
1417   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
1418   ierr = PetscFree(pa);CHKERRQ(ierr);
1419 
1420   ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
1421   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1422   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1423   for (i=0; i<nrqr; ++i) {
1424     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1425   }
1426   for (i=0; i<nrqs; ++i) {
1427     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1428   }
1429   ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr);
1430   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1431   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1432   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
1433   if (!ijonly) {
1434     for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);}
1435     ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1436     ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1437     ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
1438   }
1439 
1440 #if defined(PETSC_USE_CTABLE)
1441   for (i=0; i<ismax; i++) {
1442     ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr);
1443   }
1444 #endif
1445   ierr = PetscFree(rmap);CHKERRQ(ierr);
1446 
1447   for (i=0; i<ismax; i++) {
1448     if (!allcolumns[i]) {
1449 #if defined(PETSC_USE_CTABLE)
1450       ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr);
1451 #else
1452       ierr = PetscFree(cmap[i]);CHKERRQ(ierr);
1453 #endif
1454     }
1455   }
1456   ierr = PetscFree(cmap);CHKERRQ(ierr);
1457   ierr = PetscFree(lens);CHKERRQ(ierr);
1458 
1459   for (i=0; i<ismax; i++) {
1460     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1461     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1462   }
1463 
1464   c->ijonly = PETSC_FALSE; /* set back to the default */
1465   PetscFunctionReturn(0);
1466 }
1467 
1468