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