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