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