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