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