xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 5d0c19d75c660d4fec594a5399ec8d8ba29c54a8)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3d5b485abSSatish Balay /*
4d5b485abSSatish Balay    Routines to compute overlapping regions of a parallel MPI matrix
5d5b485abSSatish Balay   and to find submatrices that were shared across processors.
6d5b485abSSatish Balay */
770f55243SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h"
80a835dfdSSatish Balay #include "petscbt.h"
9d5b485abSSatish Balay 
10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat,PetscInt,IS *);
11b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char **,PetscInt*,PetscInt**);
12b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt **,PetscInt**,PetscInt*);
13b24ad042SBarry Smith EXTERN PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
14b24ad042SBarry Smith EXTERN PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
15d5b485abSSatish Balay 
164a2ae208SSatish Balay #undef __FUNCT__
174a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ"
18b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov)
19d5b485abSSatish Balay {
206849ba73SBarry Smith   PetscErrorCode ierr;
21d0f46423SBarry Smith   PetscInt       i,N=C->cmap->N, bs=C->rmap->bs;
2236f4e84dSSatish Balay   IS             *is_new;
2336f4e84dSSatish Balay 
243a40ed3dSBarry Smith   PetscFunctionBegin;
2582502324SSatish Balay   ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr);
26df36731bSBarry Smith   /* Convert the indices into block format */
2727f478b1SHong Zhang   ierr = ISCompressIndicesGeneral(N,bs,imax,is,is_new);CHKERRQ(ierr);
2829bbc08cSBarry Smith   if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
29d5b485abSSatish Balay   for (i=0; i<ov; ++i) {
3036f4e84dSSatish Balay     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr);
31d5b485abSSatish Balay   }
32ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);}
3327f478b1SHong Zhang   ierr = ISExpandIndicesGeneral(N,bs,imax,is_new,is);CHKERRQ(ierr);
34ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);}
35606d414cSSatish Balay   ierr = PetscFree(is_new);CHKERRQ(ierr);
363a40ed3dSBarry Smith   PetscFunctionReturn(0);
37d5b485abSSatish Balay }
38d5b485abSSatish Balay 
39d5b485abSSatish Balay /*
40d5b485abSSatish Balay   Sample message format:
41d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
420e9b0e7eSHong Zhang   to index sets is[1], is[5]
43d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
44d5b485abSSatish Balay   -----------
45d5b485abSSatish Balay   mesg [1] = 1 => is[1]
46d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
47d5b485abSSatish Balay   -----------
48d5b485abSSatish Balay   mesg [5] = 5  => is[5]
49d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
50d5b485abSSatish Balay   -----------
51d5b485abSSatish Balay   mesg [7]
520e9b0e7eSHong Zhang   mesg [n]  data(is[1])
53d5b485abSSatish Balay   -----------
54d5b485abSSatish Balay   mesg[n+1]
55d5b485abSSatish Balay   mesg[m]  data(is[5])
56d5b485abSSatish Balay   -----------
57d5b485abSSatish Balay 
58d5b485abSSatish Balay   Notes:
59d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
60d5b485abSSatish Balay   nrqr - no of requests recieved (which have to be or which have been processed
61d5b485abSSatish Balay */
624a2ae208SSatish Balay #undef __FUNCT__
634a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once"
64b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[])
65d5b485abSSatish Balay {
66df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
67*5d0c19d7SBarry Smith   const PetscInt **idx,*idx_i;
68*5d0c19d7SBarry Smith   PetscInt       *n,*w3,*w4,*rtable,**data,len;
696849ba73SBarry Smith   PetscErrorCode ierr;
70b24ad042SBarry Smith   PetscMPIInt    size,rank,tag1,tag2,*w2,*w1,nrqr;
71b24ad042SBarry Smith   PetscInt       Mbs,i,j,k,**rbuf,row,proc,nrqs,msz,**outdat,**ptr;
72b24ad042SBarry Smith   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2;
73b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1,*onodes2,*olengths2;
74f1af5d2fSBarry Smith   PetscBT        *table;
75d5b485abSSatish Balay   MPI_Comm       comm;
76d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
77d5b485abSSatish Balay   MPI_Status     *s_status,*recv_status;
78d5b485abSSatish Balay 
793a40ed3dSBarry Smith   PetscFunctionBegin;
807adad957SLisandro Dalcin   comm   = ((PetscObject)C)->comm;
81d5b485abSSatish Balay   size   = c->size;
82d5b485abSSatish Balay   rank   = c->rank;
83df36731bSBarry Smith   Mbs    = c->Mbs;
84d5b485abSSatish Balay 
85c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
86c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
87c7dd2924SSatish Balay 
88b24ad042SBarry Smith   len    = (imax+1)*sizeof(PetscInt*)+ (imax + Mbs)*sizeof(PetscInt);
8982502324SSatish Balay   ierr   = PetscMalloc(len,&idx);CHKERRQ(ierr);
90b24ad042SBarry Smith   n      = (PetscInt*)(idx + imax);
91d5b485abSSatish Balay   rtable = n + imax;
92d5b485abSSatish Balay 
93d5b485abSSatish Balay   for (i=0; i<imax; i++) {
94d5b485abSSatish Balay     ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr);
95b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
96d5b485abSSatish Balay   }
97d5b485abSSatish Balay 
98d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
99d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
100899cda47SBarry Smith     len = c->rangebs[i+1];
101d5b485abSSatish Balay     for (; j<len; j++) {
102d5b485abSSatish Balay       rtable[j] = i;
103d5b485abSSatish Balay     }
104d5b485abSSatish Balay   }
105d5b485abSSatish Balay 
106d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
107d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
108b24ad042SBarry Smith   ierr = PetscMalloc(size*2*sizeof(PetscInt)+size*2*sizeof(PetscMPIInt),&w1);CHKERRQ(ierr);/*  mesg size */
109d5b485abSSatish Balay   w2   = w1 + size;                /* if w2[i] marked, then a message to proc i*/
110b24ad042SBarry Smith   w3   = (PetscInt*) (w2 + size);   /* no of IS that needs to be sent to proc i */
111d5b485abSSatish Balay   w4   = w3 + size;                /* temp work space used in determining w1, w2, w3 */
112b24ad042SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscInt)+2*size*sizeof(PetscMPIInt));CHKERRQ(ierr); /* initialise work vector*/
113d5b485abSSatish Balay   for (i=0; i<imax; i++) {
114b24ad042SBarry Smith     ierr  = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
115d5b485abSSatish Balay     idx_i = idx[i];
116d5b485abSSatish Balay     len   = n[i];
117d5b485abSSatish Balay     for (j=0; j<len; j++) {
118d5b485abSSatish Balay       row  = idx_i[j];
1196b41c64dSBarry Smith       if (row < 0) {
120e005ede5SBarry Smith         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries");
1216b41c64dSBarry Smith       }
122d5b485abSSatish Balay       proc = rtable[row];
123d5b485abSSatish Balay       w4[proc]++;
124d5b485abSSatish Balay     }
125d5b485abSSatish Balay     for (j=0; j<size; j++){
126d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
127d5b485abSSatish Balay     }
128d5b485abSSatish Balay   }
129d5b485abSSatish Balay 
130d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
131d5b485abSSatish Balay   msz      = 0;              /* total mesg length (for all proc */
1320e9b0e7eSHong Zhang   w1[rank] = 0;              /* no mesg sent to itself */
133d5b485abSSatish Balay   w3[rank] = 0;
134d5b485abSSatish Balay   for (i=0; i<size; i++) {
135d5b485abSSatish Balay     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
136d5b485abSSatish Balay   }
137d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
138b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr);
139d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
140d5b485abSSatish Balay     if (w1[i]) {pa[j] = i; j++;}
141d5b485abSSatish Balay   }
142d5b485abSSatish Balay 
143d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
144d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
145d5b485abSSatish Balay     j      = pa[i];
146d5b485abSSatish Balay     w1[j] += w2[j] + 2*w3[j];
147d5b485abSSatish Balay     msz   += w1[j];
148d5b485abSSatish Balay   }
149d5b485abSSatish Balay 
150c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
151a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
152c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
153d5b485abSSatish Balay 
154c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
155c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr);
156d5b485abSSatish Balay 
157d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
158b24ad042SBarry Smith   len    = 2*size*sizeof(PetscInt*) + (size+msz)*sizeof(PetscInt);
15982502324SSatish Balay   ierr = PetscMalloc(len,&outdat);CHKERRQ(ierr);
160d5b485abSSatish Balay   ptr    = outdat + size;     /* Pointers to the data in outgoing buffers */
161b24ad042SBarry Smith   ierr   = PetscMemzero(outdat,2*size*sizeof(PetscInt*));CHKERRQ(ierr);
162b24ad042SBarry Smith   tmp    = (PetscInt*)(outdat + 2*size);
163d5b485abSSatish Balay   ctr    = tmp + msz;
164d5b485abSSatish Balay 
165d5b485abSSatish Balay   {
166b24ad042SBarry Smith     PetscInt *iptr = tmp,ict  = 0;
167d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
168d5b485abSSatish Balay       j         = pa[i];
169d5b485abSSatish Balay       iptr     +=  ict;
170d5b485abSSatish Balay       outdat[j] = iptr;
171d5b485abSSatish Balay       ict       = w1[j];
172d5b485abSSatish Balay     }
173d5b485abSSatish Balay   }
174d5b485abSSatish Balay 
175d5b485abSSatish Balay   /* Form the outgoing messages */
176d5b485abSSatish Balay   /*plug in the headers*/
177d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
178d5b485abSSatish Balay     j            = pa[i];
179d5b485abSSatish Balay     outdat[j][0] = 0;
180b24ad042SBarry Smith     ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
181d5b485abSSatish Balay     ptr[j]       = outdat[j] + 2*w3[j] + 1;
182d5b485abSSatish Balay   }
183d5b485abSSatish Balay 
184d5b485abSSatish Balay   /* Memory for doing local proc's work*/
185d5b485abSSatish Balay   {
186b24ad042SBarry Smith     PetscInt  *d_p;
187d5b485abSSatish Balay     char *t_p;
188d5b485abSSatish Balay 
189b24ad042SBarry Smith     len  = (imax)*(sizeof(PetscBT) + sizeof(PetscInt*)+ sizeof(PetscInt)) +
190b24ad042SBarry Smith       (Mbs)*imax*sizeof(PetscInt)  + (Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char) + 1;
19182502324SSatish Balay     ierr = PetscMalloc(len,&table);CHKERRQ(ierr);
192549d3d68SSatish Balay     ierr = PetscMemzero(table,len);CHKERRQ(ierr);
193b24ad042SBarry Smith     data = (PetscInt **)(table + imax);
194b24ad042SBarry Smith     isz  = (PetscInt  *)(data  + imax);
195b24ad042SBarry Smith     d_p  = (PetscInt  *)(isz   + imax);
196bbd702dbSSatish Balay     t_p  = (char *)(d_p   + Mbs*imax);
197bbd702dbSSatish Balay     for (i=0; i<imax; i++) {
198b6410449SSatish Balay       table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
199bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
200d5b485abSSatish Balay     }
201d5b485abSSatish Balay   }
202d5b485abSSatish Balay 
203d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
204d5b485abSSatish Balay   {
205b24ad042SBarry Smith     PetscInt     n_i,*data_i,isz_i,*outdat_j,ctr_j;
206f1af5d2fSBarry Smith     PetscBT table_i;
207d5b485abSSatish Balay 
208d5b485abSSatish Balay     for (i=0; i<imax; i++) {
209b24ad042SBarry Smith       ierr    = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
210d5b485abSSatish Balay       n_i     = n[i];
211d5b485abSSatish Balay       table_i = table[i];
212d5b485abSSatish Balay       idx_i   = idx[i];
213d5b485abSSatish Balay       data_i  = data[i];
214d5b485abSSatish Balay       isz_i   = isz[i];
215d5b485abSSatish Balay       for (j=0;  j<n_i; j++) {  /* parse the indices of each IS */
216d5b485abSSatish Balay         row  = idx_i[j];
217d5b485abSSatish Balay         proc = rtable[row];
218d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
219d5b485abSSatish Balay           ctr[proc]++;
220d5b485abSSatish Balay           *ptr[proc] = row;
221d5b485abSSatish Balay           ptr[proc]++;
222d5b485abSSatish Balay         }
223d5b485abSSatish Balay         else { /* Update the local table */
224f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
225d5b485abSSatish Balay         }
226d5b485abSSatish Balay       }
227d5b485abSSatish Balay       /* Update the headers for the current IS */
228d5b485abSSatish Balay       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
229d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
230d5b485abSSatish Balay           outdat_j        = outdat[j];
231d5b485abSSatish Balay           k               = ++outdat_j[0];
232d5b485abSSatish Balay           outdat_j[2*k]   = ctr_j;
233d5b485abSSatish Balay           outdat_j[2*k-1] = i;
234d5b485abSSatish Balay         }
235d5b485abSSatish Balay       }
236d5b485abSSatish Balay       isz[i] = isz_i;
237d5b485abSSatish Balay     }
238d5b485abSSatish Balay   }
239d5b485abSSatish Balay 
240d5b485abSSatish Balay   /*  Now  post the sends */
241b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
242d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
243d5b485abSSatish Balay     j    = pa[i];
244b24ad042SBarry Smith     ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr);
245d5b485abSSatish Balay   }
246d5b485abSSatish Balay 
247d5b485abSSatish Balay   /* No longer need the original indices*/
248d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
249d5b485abSSatish Balay     ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr);
250d5b485abSSatish Balay   }
251606d414cSSatish Balay   ierr = PetscFree(idx);CHKERRQ(ierr);
252d5b485abSSatish Balay 
253d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
254d5b485abSSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
255d5b485abSSatish Balay   }
256d5b485abSSatish Balay 
257d5b485abSSatish Balay   /* Do Local work*/
258df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
259d5b485abSSatish Balay 
260d5b485abSSatish Balay   /* Receive messages*/
261b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr);
2620c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);}
263d5b485abSSatish Balay 
264b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr);
2650c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);}
266d5b485abSSatish Balay 
267d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
268606d414cSSatish Balay   ierr = PetscFree(outdat);CHKERRQ(ierr);
269606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
270d5b485abSSatish Balay 
271b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&xdata);CHKERRQ(ierr);
272b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt),&isz1);CHKERRQ(ierr);
273df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
274606d414cSSatish Balay   ierr = PetscFree(rbuf);CHKERRQ(ierr);
275d5b485abSSatish Balay 
276d5b485abSSatish Balay   /* Send the data back*/
277d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
278d5b485abSSatish Balay   {
279b24ad042SBarry Smith     PetscMPIInt *rw1;
280d5b485abSSatish Balay 
281b24ad042SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&rw1);CHKERRQ(ierr);
282b24ad042SBarry Smith     ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr);
283c7dd2924SSatish Balay 
284d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
285d5b485abSSatish Balay       proc      = recv_status[i].MPI_SOURCE;
286e005ede5SBarry Smith       if (proc != onodes1[i]) SETERRQ(PETSC_ERR_PLIB,"MPI_SOURCE mismatch");
287d5b485abSSatish Balay       rw1[proc] = isz1[i];
288d5b485abSSatish Balay     }
289d5b485abSSatish Balay 
290c7dd2924SSatish Balay     ierr = PetscFree(onodes1);CHKERRQ(ierr);
291c7dd2924SSatish Balay     ierr = PetscFree(olengths1);CHKERRQ(ierr);
292d5b485abSSatish Balay 
293c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
294c7dd2924SSatish Balay     ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr);
29503512dc5SSatish Balay     ierr = PetscFree(rw1);CHKERRQ(ierr);
296c7dd2924SSatish Balay   }
297c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
298c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr);
299d5b485abSSatish Balay 
300d5b485abSSatish Balay   /*  Now  post the sends */
301b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
302d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
303d5b485abSSatish Balay     j    = recv_status[i].MPI_SOURCE;
304b24ad042SBarry Smith     ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr);
305d5b485abSSatish Balay   }
306d5b485abSSatish Balay 
307d5b485abSSatish Balay   /* receive work done on other processors*/
308d5b485abSSatish Balay   {
309b24ad042SBarry Smith     PetscMPIInt idex;
310b24ad042SBarry Smith     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
311f1af5d2fSBarry Smith     PetscBT     table_i;
312d5b485abSSatish Balay     MPI_Status  *status2;
313d5b485abSSatish Balay 
314169449a0SSatish Balay     ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr);
315d5b485abSSatish Balay 
316d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
317999d9058SBarry Smith       ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr);
318d5b485abSSatish Balay       /* Process the message*/
319999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
320d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
321999d9058SBarry Smith       jmax    = rbuf2[idex][0];
322d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
323d5b485abSSatish Balay         max     = rbuf2_i[2*j];
324d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
325d5b485abSSatish Balay         isz_i   = isz[is_no];
326d5b485abSSatish Balay         data_i  = data[is_no];
327d5b485abSSatish Balay         table_i = table[is_no];
328d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
329d5b485abSSatish Balay           row = rbuf2_i[ct1];
330f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
331d5b485abSSatish Balay         }
332d5b485abSSatish Balay         isz[is_no] = isz_i;
333d5b485abSSatish Balay       }
334d5b485abSSatish Balay     }
3350c468ba9SBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);}
336606d414cSSatish Balay     ierr = PetscFree(status2);CHKERRQ(ierr);
337d5b485abSSatish Balay   }
338d5b485abSSatish Balay 
339d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
340029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr);
341d5b485abSSatish Balay   }
342d5b485abSSatish Balay 
343c7dd2924SSatish Balay 
344c7dd2924SSatish Balay   ierr = PetscFree(onodes2);CHKERRQ(ierr);
345c7dd2924SSatish Balay   ierr = PetscFree(olengths2);CHKERRQ(ierr);
346c7dd2924SSatish Balay 
347606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
348606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
349606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
350606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
351606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
352606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
353606d414cSSatish Balay   ierr = PetscFree(table);CHKERRQ(ierr);
354606d414cSSatish Balay   ierr = PetscFree(s_status);CHKERRQ(ierr);
355606d414cSSatish Balay   ierr = PetscFree(recv_status);CHKERRQ(ierr);
356606d414cSSatish Balay   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
357606d414cSSatish Balay   ierr = PetscFree(xdata);CHKERRQ(ierr);
358606d414cSSatish Balay   ierr = PetscFree(isz1);CHKERRQ(ierr);
3593a40ed3dSBarry Smith   PetscFunctionReturn(0);
360d5b485abSSatish Balay }
361d5b485abSSatish Balay 
3624a2ae208SSatish Balay #undef __FUNCT__
3634a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local"
364d5b485abSSatish Balay /*
365df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
366d5b485abSSatish Balay        the work on the local processor.
367d5b485abSSatish Balay 
368d5b485abSSatish Balay      Inputs:
369df36731bSBarry Smith       C      - MAT_MPIBAIJ;
370d5b485abSSatish Balay       imax - total no of index sets processed at a time;
371df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
372d5b485abSSatish Balay 
373d5b485abSSatish Balay      Output:
3740e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
375d5b485abSSatish Balay                to each index set;
376d5b485abSSatish Balay       data   - pointer to the solutions
377d5b485abSSatish Balay */
378b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data)
379d5b485abSSatish Balay {
380df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
381d5b485abSSatish Balay   Mat         A = c->A,B = c->B;
382df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
383b24ad042SBarry Smith   PetscInt    start,end,val,max,rstart,cstart,*ai,*aj;
384b24ad042SBarry Smith   PetscInt    *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
385f1af5d2fSBarry Smith   PetscBT     table_i;
386d5b485abSSatish Balay 
3873a40ed3dSBarry Smith   PetscFunctionBegin;
388899cda47SBarry Smith   rstart = c->rstartbs;
389899cda47SBarry Smith   cstart = c->cstartbs;
390d5b485abSSatish Balay   ai     = a->i;
391df36731bSBarry Smith   aj     = a->j;
392d5b485abSSatish Balay   bi     = b->i;
393df36731bSBarry Smith   bj     = b->j;
394d5b485abSSatish Balay   garray = c->garray;
395d5b485abSSatish Balay 
396d5b485abSSatish Balay 
397d5b485abSSatish Balay   for (i=0; i<imax; i++) {
398d5b485abSSatish Balay     data_i  = data[i];
399d5b485abSSatish Balay     table_i = table[i];
400d5b485abSSatish Balay     isz_i   = isz[i];
401d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
402d5b485abSSatish Balay       row   = data_i[j] - rstart;
403d5b485abSSatish Balay       start = ai[row];
404d5b485abSSatish Balay       end   = ai[row+1];
405d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
406df36731bSBarry Smith         val = aj[k] + cstart;
407f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
408d5b485abSSatish Balay       }
409d5b485abSSatish Balay       start = bi[row];
410d5b485abSSatish Balay       end   = bi[row+1];
411d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
412df36731bSBarry Smith         val = garray[bj[k]];
413f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
414d5b485abSSatish Balay       }
415d5b485abSSatish Balay     }
416d5b485abSSatish Balay     isz[i] = isz_i;
417d5b485abSSatish Balay   }
4183a40ed3dSBarry Smith   PetscFunctionReturn(0);
419d5b485abSSatish Balay }
4204a2ae208SSatish Balay #undef __FUNCT__
4214a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive"
422d5b485abSSatish Balay /*
423df36731bSBarry Smith       MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages,
424d5b485abSSatish Balay          and return the output
425d5b485abSSatish Balay 
426d5b485abSSatish Balay          Input:
427d5b485abSSatish Balay            C    - the matrix
428d5b485abSSatish Balay            nrqr - no of messages being processed.
429d5b485abSSatish Balay            rbuf - an array of pointers to the recieved requests
430d5b485abSSatish Balay 
431d5b485abSSatish Balay          Output:
432d5b485abSSatish Balay            xdata - array of messages to be sent back
433d5b485abSSatish Balay            isz1  - size of each message
434d5b485abSSatish Balay 
435a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
436d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
4370e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of
438d5b485abSSatish Balay memory is used.
439d5b485abSSatish Balay 
440d5b485abSSatish Balay */
441b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
442d5b485abSSatish Balay {
443df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
444d5b485abSSatish Balay   Mat            A = c->A,B = c->B;
445df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
4466849ba73SBarry Smith   PetscErrorCode ierr;
447b24ad042SBarry Smith   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
448b24ad042SBarry Smith   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
449b24ad042SBarry Smith   PetscInt       val,max1,max2,rank,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
450b24ad042SBarry Smith   PetscInt       *rbuf_i,kmax,rbuf_0;
451f1af5d2fSBarry Smith   PetscBT        xtable;
452d5b485abSSatish Balay 
4533a40ed3dSBarry Smith   PetscFunctionBegin;
454d5b485abSSatish Balay   rank   = c->rank;
455df36731bSBarry Smith   Mbs    = c->Mbs;
456899cda47SBarry Smith   rstart = c->rstartbs;
457899cda47SBarry Smith   cstart = c->cstartbs;
458d5b485abSSatish Balay   ai     = a->i;
459df36731bSBarry Smith   aj     = a->j;
460d5b485abSSatish Balay   bi     = b->i;
461df36731bSBarry Smith   bj     = b->j;
462d5b485abSSatish Balay   garray = c->garray;
463d5b485abSSatish Balay 
464d5b485abSSatish Balay 
465d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
466d5b485abSSatish Balay     rbuf_i  =  rbuf[i];
467d5b485abSSatish Balay     rbuf_0  =  rbuf_i[0];
468d5b485abSSatish Balay     ct     += rbuf_0;
469d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; }
470d5b485abSSatish Balay   }
471d5b485abSSatish Balay 
472701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
473701b8814SBarry Smith   else        max1 = 1;
474d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
475b24ad042SBarry Smith   ierr         = PetscMalloc(mem_estimate*sizeof(PetscInt),&xdata[0]);CHKERRQ(ierr);
476d5b485abSSatish Balay   ++no_malloc;
4776831982aSBarry Smith   ierr         = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr);
478b24ad042SBarry Smith   ierr         = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr);
479d5b485abSSatish Balay 
480d5b485abSSatish Balay   ct3 = 0;
481d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
482d5b485abSSatish Balay     rbuf_i =  rbuf[i];
483d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
484d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
485d5b485abSSatish Balay     ct2    =  ct1;
486d5b485abSSatish Balay     ct3    += ct1;
487d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
4886831982aSBarry Smith       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
489d5b485abSSatish Balay       oct2 = ct2;
490d5b485abSSatish Balay       kmax = rbuf_i[2*j];
491d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
492d5b485abSSatish Balay         row = rbuf_i[ct1];
493f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
494d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
495b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
496b24ad042SBarry Smith             ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
497b24ad042SBarry Smith             ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
498606d414cSSatish Balay             ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
499d5b485abSSatish Balay             xdata[0]     = tmp;
500d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
501d5b485abSSatish Balay             for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
502d5b485abSSatish Balay           }
503d5b485abSSatish Balay           xdata[i][ct2++] = row;
504d5b485abSSatish Balay           ct3++;
505d5b485abSSatish Balay         }
506d5b485abSSatish Balay       }
507d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
508d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
509d5b485abSSatish Balay         start = ai[row];
510d5b485abSSatish Balay         end   = ai[row+1];
511d5b485abSSatish Balay         for (l=start; l<end; l++) {
512df36731bSBarry Smith           val = aj[l] + cstart;
513f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
514d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
515b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
516b24ad042SBarry Smith               ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
517b24ad042SBarry Smith               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
518606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
519d5b485abSSatish Balay               xdata[0]     = tmp;
520d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
521d5b485abSSatish Balay               for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
522d5b485abSSatish Balay             }
523d5b485abSSatish Balay             xdata[i][ct2++] = val;
524d5b485abSSatish Balay             ct3++;
525d5b485abSSatish Balay           }
526d5b485abSSatish Balay         }
527d5b485abSSatish Balay         start = bi[row];
528d5b485abSSatish Balay         end   = bi[row+1];
529d5b485abSSatish Balay         for (l=start; l<end; l++) {
530df36731bSBarry Smith           val = garray[bj[l]];
531f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
532d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
533b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
534b24ad042SBarry Smith               ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
535b24ad042SBarry Smith               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
536606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
537d5b485abSSatish Balay               xdata[0]     = tmp;
538d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
539d5b485abSSatish Balay               for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
540d5b485abSSatish Balay             }
541d5b485abSSatish Balay             xdata[i][ct2++] = val;
542d5b485abSSatish Balay             ct3++;
543d5b485abSSatish Balay           }
544d5b485abSSatish Balay         }
545d5b485abSSatish Balay       }
546d5b485abSSatish Balay       /* Update the header*/
547d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
548d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
549d5b485abSSatish Balay     }
550d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
551d5b485abSSatish Balay     xdata[i+1]  = xdata[i] + ct2;
552d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
553d5b485abSSatish Balay   }
5546831982aSBarry Smith   ierr = PetscBTDestroy(xtable);CHKERRQ(ierr);
5551e2582c4SBarry Smith   ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr);
5563a40ed3dSBarry Smith   PetscFunctionReturn(0);
557d5b485abSSatish Balay }
558d5b485abSSatish Balay 
559b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *);
560a2feddc5SSatish Balay 
5614a2ae208SSatish Balay #undef __FUNCT__
5624a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ"
563b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
564d5b485abSSatish Balay {
56536f4e84dSSatish Balay   IS             *isrow_new,*iscol_new;
566cf4f527aSSatish Balay   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
5676849ba73SBarry Smith   PetscErrorCode ierr;
568d0f46423SBarry Smith   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,N=C->cmap->N,bs=C->rmap->bs;
569a2feddc5SSatish Balay 
5703a40ed3dSBarry Smith   PetscFunctionBegin;
5713a40ed3dSBarry Smith   /* The compression and expansion should be avoided. Does'nt point
572a2feddc5SSatish Balay      out errors might change the indices hence buggey */
573a2feddc5SSatish Balay 
574b0a32e0cSBarry Smith   ierr = PetscMalloc(2*(ismax+1)*sizeof(IS),&isrow_new);CHKERRQ(ierr);
57536f4e84dSSatish Balay   iscol_new = isrow_new + ismax;
57627f478b1SHong Zhang   ierr = ISCompressIndicesSorted(N,bs,ismax,isrow,isrow_new);CHKERRQ(ierr);
57727f478b1SHong Zhang   ierr = ISCompressIndicesSorted(N,bs,ismax,iscol,iscol_new);CHKERRQ(ierr);
578cf4f527aSSatish Balay 
579cf4f527aSSatish Balay   /* Allocate memory to hold all the submatrices */
580cf4f527aSSatish Balay   if (scall != MAT_REUSE_MATRIX) {
58182502324SSatish Balay     ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr);
582cf4f527aSSatish Balay   }
583cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
584b24ad042SBarry Smith   nmax          = 20*1000000 / (c->Nbs * sizeof(PetscInt));
585329f5518SBarry Smith   if (!nmax) nmax = 1;
586cf4f527aSSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
587cf4f527aSSatish Balay 
588653e4784SBarry Smith   /* Make sure every processor loops through the nstages */
5897adad957SLisandro Dalcin   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);CHKERRQ(ierr);
590cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
591cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
592cf4f527aSSatish Balay     else if (pos == ismax) max_no = 0;
593cf4f527aSSatish Balay     else                   max_no = ismax-pos;
594cf4f527aSSatish Balay     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr);
595cf4f527aSSatish Balay     pos += max_no;
596cf4f527aSSatish Balay   }
59736f4e84dSSatish Balay 
59836f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
599ca161407SBarry Smith     ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr);
600ca161407SBarry Smith     ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr);
60136f4e84dSSatish Balay   }
602606d414cSSatish Balay   ierr = PetscFree(isrow_new);CHKERRQ(ierr);
6033a40ed3dSBarry Smith   PetscFunctionReturn(0);
604a2feddc5SSatish Balay }
605a2feddc5SSatish Balay 
606233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
6074a2ae208SSatish Balay #undef __FUNCT__
6084a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc"
609e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
610233c2ff6SSatish Balay {
611e005ede5SBarry Smith   PetscInt    nGlobalNd = proc_gnode[size];
612b9f7ace7SBarry Smith   PetscMPIInt fproc = PetscMPIIntCast( (PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)));
613233c2ff6SSatish Balay 
614233c2ff6SSatish Balay   PetscFunctionBegin;
61523ce1328SBarry Smith   if (fproc > size) fproc = size;
616e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
617e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
618233c2ff6SSatish Balay     else                         fproc++;
619233c2ff6SSatish Balay   }
620e005ede5SBarry Smith   *rank = fproc;
621233c2ff6SSatish Balay   PetscFunctionReturn(0);
622233c2ff6SSatish Balay }
623233c2ff6SSatish Balay #endif
624233c2ff6SSatish Balay 
625a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
626b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
6274a2ae208SSatish Balay #undef __FUNCT__
6284a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local"
629b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
630a2feddc5SSatish Balay {
631df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
632cf4f527aSSatish Balay   Mat            A = c->A;
633df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
634*5d0c19d7SBarry Smith   const PetscInt **irow,**icol,*irow_i;
635*5d0c19d7SBarry Smith   PetscInt       *nrow,*ncol,*w3,*w4,start;
6366849ba73SBarry Smith   PetscErrorCode ierr;
6379405f653SBarry Smith   PetscMPIInt    size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc;
6389405f653SBarry Smith   PetscInt       **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row;
639b24ad042SBarry Smith   PetscInt       nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol;
640b24ad042SBarry Smith   PetscInt       **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
641*5d0c19d7SBarry Smith   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
642b24ad042SBarry Smith   PetscInt       len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
643d0f46423SBarry Smith   PetscInt       bs=C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
644899cda47SBarry Smith   PetscInt       cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
645899cda47SBarry Smith   PetscInt       *bmap = c->garray,ctmp,rstart=c->rstartbs;
646d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
647d5b485abSSatish Balay   MPI_Request    *r_waits4,*s_waits3,*s_waits4;
648d5b485abSSatish Balay   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
649d5b485abSSatish Balay   MPI_Status     *r_status3,*r_status4,*s_status4;
650d5b485abSSatish Balay   MPI_Comm       comm;
6513eda8832SBarry Smith   MatScalar      **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB;
6523eda8832SBarry Smith   MatScalar      *a_a=a->a,*b_a=b->a;
6530f5bd95cSBarry Smith   PetscTruth     flag;
654b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1;
655c7dd2924SSatish Balay 
65680d608c0SSatish Balay #if defined (PETSC_USE_CTABLE)
657b24ad042SBarry Smith   PetscInt tt;
658233c2ff6SSatish Balay   PetscTable  *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1;
659233c2ff6SSatish Balay #else
660b24ad042SBarry Smith   PetscInt         **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
661233c2ff6SSatish Balay #endif
662d5b485abSSatish Balay 
6633a40ed3dSBarry Smith   PetscFunctionBegin;
6647adad957SLisandro Dalcin   comm   = ((PetscObject)C)->comm;
6657adad957SLisandro Dalcin   tag0   = ((PetscObject)C)->tag;
666d5b485abSSatish Balay   size   = c->size;
667d5b485abSSatish Balay   rank   = c->rank;
668d5b485abSSatish Balay 
66934e6ae18SSatish Balay   /* Get some new tags to keep the communication clean */
67034e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
67134e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
67234e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
67334e6ae18SSatish Balay 
674d5b485abSSatish Balay   /* Check if the col indices are sorted */
675d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
676888f2ed8SSatish Balay     ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr);
67729bbc08cSBarry Smith     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
678d5b485abSSatish Balay   }
679d5b485abSSatish Balay 
680b24ad042SBarry Smith   len    = (2*ismax+1)*(sizeof(PetscInt*)+ sizeof(PetscInt));
681233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE)
682b24ad042SBarry Smith   len    += (Mbs+1)*sizeof(PetscInt);
683233c2ff6SSatish Balay #endif
68482502324SSatish Balay   ierr = PetscMalloc(len,&irow);CHKERRQ(ierr);
685d5b485abSSatish Balay   icol = irow + ismax;
686b24ad042SBarry Smith   nrow = (PetscInt*)(icol + ismax);
687d5b485abSSatish Balay   ncol = nrow + ismax;
688233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE)
689d5b485abSSatish Balay   rtable = ncol + ismax;
690d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
691d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
692d5b485abSSatish Balay     jmax = c->rowners[i+1];
693d5b485abSSatish Balay     for (; j<jmax; j++) {
694d5b485abSSatish Balay       rtable[j] = i;
695d5b485abSSatish Balay     }
696d5b485abSSatish Balay   }
697233c2ff6SSatish Balay #endif
698233c2ff6SSatish Balay 
699233c2ff6SSatish Balay   for (i=0; i<ismax; i++) {
700233c2ff6SSatish Balay     ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
701233c2ff6SSatish Balay     ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
702233c2ff6SSatish Balay     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
703233c2ff6SSatish Balay     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
704233c2ff6SSatish Balay   }
705d5b485abSSatish Balay 
706d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg,and buffer space
707d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
708b24ad042SBarry Smith   ierr = PetscMalloc(size*2*sizeof(PetscInt) + size*2*sizeof(PetscMPIInt),&w1);CHKERRQ(ierr); /* mesg size */
709d5b485abSSatish Balay   w2   = w1 + size;               /* if w2[i] marked, then a message to proc i*/
710b24ad042SBarry Smith   w3   = (PetscInt*)(w2 + size);  /* no of IS that needs to be sent to proc i */
711d5b485abSSatish Balay   w4   = w3 + size;                /* temp work space used in determining w1, w2, w3 */
712b24ad042SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscInt)+2*size*sizeof(PetscMPIInt));CHKERRQ(ierr); /* initialise work vector*/
713d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
714b24ad042SBarry Smith     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
715d5b485abSSatish Balay     jmax   = nrow[i];
716d5b485abSSatish Balay     irow_i = irow[i];
717d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
718d5b485abSSatish Balay       row  = irow_i[j];
719233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
720899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
721233c2ff6SSatish Balay #else
722d5b485abSSatish Balay       proc = rtable[row];
723233c2ff6SSatish Balay #endif
724d5b485abSSatish Balay       w4[proc]++;
725d5b485abSSatish Balay     }
726d5b485abSSatish Balay     for (j=0; j<size; j++) {
727d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
728d5b485abSSatish Balay     }
729d5b485abSSatish Balay   }
730d5b485abSSatish Balay 
731d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
732e757655aSSatish Balay   msz      = 0;              /* total mesg length for all proc */
733d5b485abSSatish Balay   w1[rank] = 0;              /* no mesg sent to intself */
734d5b485abSSatish Balay   w3[rank] = 0;
735d5b485abSSatish Balay   for (i=0; i<size; i++) {
736d5b485abSSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
737d5b485abSSatish Balay   }
738b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/
739d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
740d5b485abSSatish Balay     if (w1[i]) { pa[j] = i; j++; }
741d5b485abSSatish Balay   }
742d5b485abSSatish Balay 
743d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
744d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
745d5b485abSSatish Balay     j     = pa[i];
746d5b485abSSatish Balay     w1[j] += w2[j] + 2* w3[j];
747d5b485abSSatish Balay     msz   += w1[j];
748d5b485abSSatish Balay   }
749d5b485abSSatish Balay 
750c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
751a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
752c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
753d5b485abSSatish Balay 
754c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
755c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
756c7dd2924SSatish Balay 
757c7dd2924SSatish Balay   ierr = PetscFree(onodes1);CHKERRQ(ierr);
758c7dd2924SSatish Balay   ierr = PetscFree(olengths1);CHKERRQ(ierr);
759d5b485abSSatish Balay 
760d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
761b24ad042SBarry Smith   len  = 2*size*sizeof(PetscInt*) + 2*msz*sizeof(PetscInt) + size*sizeof(PetscInt);
76282502324SSatish Balay   ierr = PetscMalloc(len,&sbuf1);CHKERRQ(ierr);
763d5b485abSSatish Balay   ptr  = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
764b24ad042SBarry Smith   ierr = PetscMemzero(sbuf1,2*size*sizeof(PetscInt*));CHKERRQ(ierr);
765d5b485abSSatish Balay   /* allocate memory for outgoing data + buf to receive the first reply */
766b24ad042SBarry Smith   tmp  = (PetscInt*)(ptr + size);
767d5b485abSSatish Balay   ctr  = tmp + 2*msz;
768d5b485abSSatish Balay 
769d5b485abSSatish Balay   {
770b24ad042SBarry Smith     PetscInt *iptr = tmp,ict = 0;
771d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
772d5b485abSSatish Balay       j         = pa[i];
773d5b485abSSatish Balay       iptr     += ict;
774d5b485abSSatish Balay       sbuf1[j]  = iptr;
775d5b485abSSatish Balay       ict       = w1[j];
776d5b485abSSatish Balay     }
777d5b485abSSatish Balay   }
778d5b485abSSatish Balay 
779d5b485abSSatish Balay   /* Form the outgoing messages */
780d5b485abSSatish Balay   /* Initialise the header space */
781d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
782d5b485abSSatish Balay     j           = pa[i];
783d5b485abSSatish Balay     sbuf1[j][0] = 0;
784b24ad042SBarry Smith     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
785d5b485abSSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
786d5b485abSSatish Balay   }
787d5b485abSSatish Balay 
788d5b485abSSatish Balay   /* Parse the isrow and copy data into outbuf */
789d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
790b24ad042SBarry Smith     ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
791d5b485abSSatish Balay     irow_i = irow[i];
792d5b485abSSatish Balay     jmax   = nrow[i];
793d5b485abSSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
794d5b485abSSatish Balay       row  = irow_i[j];
795233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
796899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
797233c2ff6SSatish Balay #else
798d5b485abSSatish Balay       proc = rtable[row];
799233c2ff6SSatish Balay #endif
800d5b485abSSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
801d5b485abSSatish Balay         ctr[proc]++;
802d5b485abSSatish Balay         *ptr[proc] = row;
803d5b485abSSatish Balay         ptr[proc]++;
804d5b485abSSatish Balay       }
805d5b485abSSatish Balay     }
806d5b485abSSatish Balay     /* Update the headers for the current IS */
807d5b485abSSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
80806ef35abSSatish Balay       if ((ctr_j = ctr[j])) {
809d5b485abSSatish Balay         sbuf1_j        = sbuf1[j];
810d5b485abSSatish Balay         k              = ++sbuf1_j[0];
811d5b485abSSatish Balay         sbuf1_j[2*k]   = ctr_j;
812d5b485abSSatish Balay         sbuf1_j[2*k-1] = i;
813d5b485abSSatish Balay       }
814d5b485abSSatish Balay     }
815d5b485abSSatish Balay   }
816d5b485abSSatish Balay 
817d5b485abSSatish Balay   /*  Now  post the sends */
818b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
819d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
820d5b485abSSatish Balay     j = pa[i];
821b24ad042SBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
822d5b485abSSatish Balay   }
823d5b485abSSatish Balay 
824d5b485abSSatish Balay   /* Post Recieves to capture the buffer size */
825b0a32e0cSBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr);
826b24ad042SBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr);
827d5b485abSSatish Balay   rbuf2[0] = tmp + msz;
828d5b485abSSatish Balay   for (i=1; i<nrqs; ++i) {
829d5b485abSSatish Balay     j        = pa[i];
830d5b485abSSatish Balay     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
831d5b485abSSatish Balay   }
832d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
833d5b485abSSatish Balay     j    = pa[i];
834b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
835d5b485abSSatish Balay   }
836d5b485abSSatish Balay 
837d5b485abSSatish Balay   /* Send to other procs the buf size they should allocate */
838d5b485abSSatish Balay 
839d5b485abSSatish Balay   /* Receive messages*/
840b0a32e0cSBarry Smith   ierr       = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
841b0a32e0cSBarry Smith   ierr       = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr);
842b24ad042SBarry Smith   len        = 2*nrqr*sizeof(PetscInt) + (nrqr+1)*sizeof(PetscInt*);
84382502324SSatish Balay   ierr       = PetscMalloc(len,&sbuf2);CHKERRQ(ierr);
844b24ad042SBarry Smith   req_size   = (PetscInt*)(sbuf2 + nrqr);
845d5b485abSSatish Balay   req_source = req_size + nrqr;
846d5b485abSSatish Balay 
847d5b485abSSatish Balay   {
848df36731bSBarry Smith     Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
849b24ad042SBarry Smith     PetscInt    *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
850d5b485abSSatish Balay 
851d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
852999d9058SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
853999d9058SBarry Smith       req_size[idex] = 0;
854999d9058SBarry Smith       rbuf1_i         = rbuf1[idex];
855d5b485abSSatish Balay       start           = 2*rbuf1_i[0] + 1;
856b24ad042SBarry Smith       ierr            = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
857b24ad042SBarry Smith       ierr            = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr);
858999d9058SBarry Smith       sbuf2_i         = sbuf2[idex];
859d5b485abSSatish Balay       for (j=start; j<end; j++) {
860d5b485abSSatish Balay         id               = rbuf1_i[j] - rstart;
861d5b485abSSatish Balay         ncols            = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
862d5b485abSSatish Balay         sbuf2_i[j]       = ncols;
863999d9058SBarry Smith         req_size[idex] += ncols;
864d5b485abSSatish Balay       }
865999d9058SBarry Smith       req_source[idex] = r_status1[i].MPI_SOURCE;
866d5b485abSSatish Balay       /* form the header */
867999d9058SBarry Smith       sbuf2_i[0]   = req_size[idex];
868d5b485abSSatish Balay       for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
869b24ad042SBarry Smith       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
870d5b485abSSatish Balay     }
871d5b485abSSatish Balay   }
872606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
873606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
874d5b485abSSatish Balay 
875d5b485abSSatish Balay   /*  recv buffer sizes */
876d5b485abSSatish Balay   /* Receive messages*/
877d5b485abSSatish Balay 
878b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr);
879b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr);
880b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr);
881b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr);
882b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr);
883d5b485abSSatish Balay 
884d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
885999d9058SBarry Smith     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
886b24ad042SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr);
887999d9058SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr);
888b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
889b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
890d5b485abSSatish Balay   }
891606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
892606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
893d5b485abSSatish Balay 
894d5b485abSSatish Balay   /* Wait on sends1 and sends2 */
895b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr);
896b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr);
897d5b485abSSatish Balay 
8980c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
8990c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
900606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
901606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
902606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
903606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
904d5b485abSSatish Balay 
905d5b485abSSatish Balay   /* Now allocate buffers for a->j, and send them off */
906b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr);
907d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
908b24ad042SBarry Smith   ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr);
909d5b485abSSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
910d5b485abSSatish Balay 
911b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr);
912d5b485abSSatish Balay   {
913d5b485abSSatish Balay      for (i=0; i<nrqr; i++) {
914d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
915d5b485abSSatish Balay       sbuf_aj_i = sbuf_aj[i];
916d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0] + 1;
917d5b485abSSatish Balay       ct2       = 0;
918d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
919d5b485abSSatish Balay         kmax = rbuf1[i][2*j];
920d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
921e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
922d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
923d5b485abSSatish Balay           ncols  = nzA + nzB;
924d5b485abSSatish Balay           cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
925d5b485abSSatish Balay 
926d5b485abSSatish Balay           /* load the column indices for this row into cols*/
927d5b485abSSatish Balay           cols  = sbuf_aj_i + ct2;
928d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
929d5b485abSSatish Balay             if ((ctmp = bmap[cworkB[l]]) < cstart)  cols[l] = ctmp;
930d5b485abSSatish Balay             else break;
931d5b485abSSatish Balay           }
932d5b485abSSatish Balay           imark = l;
933d5b485abSSatish Balay           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
934d5b485abSSatish Balay           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
935d5b485abSSatish Balay           ct2 += ncols;
936d5b485abSSatish Balay         }
937d5b485abSSatish Balay       }
938b24ad042SBarry Smith       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
939d5b485abSSatish Balay     }
940d5b485abSSatish Balay   }
941b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr);
942b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr);
943d5b485abSSatish Balay 
944d5b485abSSatish Balay   /* Allocate buffers for a->a, and send them off */
94582502324SSatish Balay   ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr);
946d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
94782502324SSatish Balay   ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr);
948a2feddc5SSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
949d5b485abSSatish Balay 
950b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr);
951d5b485abSSatish Balay   {
952d5b485abSSatish Balay     for (i=0; i<nrqr; i++) {
953d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
954d5b485abSSatish Balay       sbuf_aa_i = sbuf_aa[i];
955d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0]+1;
956d5b485abSSatish Balay       ct2       = 0;
957d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
958d5b485abSSatish Balay         kmax = rbuf1_i[2*j];
959d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
960e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
961d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
962d5b485abSSatish Balay           ncols  = nzA + nzB;
963d5b485abSSatish Balay           cworkA = a_j + a_i[row];     cworkB = b_j + b_i[row];
964a2feddc5SSatish Balay           vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2;
965d5b485abSSatish Balay 
966d5b485abSSatish Balay           /* load the column values for this row into vals*/
9675b83ace0SSatish Balay           vals  = sbuf_aa_i+ct2*bs2;
968d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
9693a40ed3dSBarry Smith             if ((bmap[cworkB[l]]) < cstart) {
9703eda8832SBarry Smith               ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9713a40ed3dSBarry Smith             }
972d5b485abSSatish Balay             else break;
973d5b485abSSatish Balay           }
974d5b485abSSatish Balay           imark = l;
9753a40ed3dSBarry Smith           for (l=0; l<nzA; l++) {
9763eda8832SBarry Smith             ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9773a40ed3dSBarry Smith           }
9783a40ed3dSBarry Smith           for (l=imark; l<nzB; l++) {
9793eda8832SBarry Smith             ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9803a40ed3dSBarry Smith           }
981d5b485abSSatish Balay           ct2 += ncols;
982d5b485abSSatish Balay         }
983d5b485abSSatish Balay       }
9843eda8832SBarry Smith       ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
985d5b485abSSatish Balay     }
986d5b485abSSatish Balay   }
987b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr);
988b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr);
989606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
990d5b485abSSatish Balay 
991d5b485abSSatish Balay   /* Form the matrix */
992d5b485abSSatish Balay   /* create col map */
993d5b485abSSatish Balay   {
994*5d0c19d7SBarry Smith     const PetscInt *icol_i;
995233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
996233c2ff6SSatish Balay     /* Create row map*/
997b0a32e0cSBarry Smith     ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&colmaps);CHKERRQ(ierr);
998ff0794f7SSatish Balay     for (i=0; i<ismax; i++) {
999ff0794f7SSatish Balay       ierr = PetscTableCreate(ncol[i]+1,&colmaps[i]);CHKERRQ(ierr);
1000233c2ff6SSatish Balay     }
1001233c2ff6SSatish Balay #else
1002b24ad042SBarry Smith     len     = (1+ismax)*sizeof(PetscInt*)+ ismax*c->Nbs*sizeof(PetscInt);
100382502324SSatish Balay     ierr    = PetscMalloc(len,&cmap);CHKERRQ(ierr);
1004b24ad042SBarry Smith     cmap[0] = (PetscInt*)(cmap + ismax);
1005b24ad042SBarry Smith     ierr    = PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(PetscInt));CHKERRQ(ierr);
1006a2feddc5SSatish Balay     for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; }
1007233c2ff6SSatish Balay #endif
1008d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1009d5b485abSSatish Balay       jmax   = ncol[i];
1010d5b485abSSatish Balay       icol_i = icol[i];
1011233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1012233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
1013233c2ff6SSatish Balay       for (j=0; j<jmax; j++) {
1014001aedefSBarry Smith 	ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr);
1015233c2ff6SSatish Balay       }
1016233c2ff6SSatish Balay #else
1017d5b485abSSatish Balay       cmap_i = cmap[i];
1018d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1019d5b485abSSatish Balay         cmap_i[icol_i[j]] = j+1;
1020d5b485abSSatish Balay       }
1021233c2ff6SSatish Balay #endif
1022d5b485abSSatish Balay     }
1023d5b485abSSatish Balay   }
1024d5b485abSSatish Balay 
1025d5b485abSSatish Balay   /* Create lens which is required for MatCreate... */
1026d5b485abSSatish Balay   for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
1027b24ad042SBarry Smith   len     = (1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt);
102882502324SSatish Balay   ierr    = PetscMalloc(len,&lens);CHKERRQ(ierr);
1029b24ad042SBarry Smith   lens[0] = (PetscInt*)(lens + ismax);
1030b24ad042SBarry Smith   ierr    = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr);
1031d5b485abSSatish Balay   for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
1032d5b485abSSatish Balay 
1033d5b485abSSatish Balay   /* Update lens from local data */
1034d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1035d5b485abSSatish Balay     jmax   = nrow[i];
1036233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1037233c2ff6SSatish Balay     lcol1_gcol1 = colmaps[i];
1038233c2ff6SSatish Balay #else
1039d5b485abSSatish Balay     cmap_i = cmap[i];
1040233c2ff6SSatish Balay #endif
1041d5b485abSSatish Balay     irow_i = irow[i];
1042d5b485abSSatish Balay     lens_i = lens[i];
1043d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1044d5b485abSSatish Balay       row  = irow_i[j];
1045233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1046899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1047233c2ff6SSatish Balay #else
1048d5b485abSSatish Balay       proc = rtable[row];
1049233c2ff6SSatish Balay #endif
1050d5b485abSSatish Balay       if (proc == rank) {
105136f4e84dSSatish Balay         /* Get indices from matA and then from matB */
105236f4e84dSSatish Balay         row    = row - rstart;
105336f4e84dSSatish Balay         nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
105436f4e84dSSatish Balay         cworkA =  a_j + a_i[row]; cworkB = b_j + b_i[row];
1055233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1056233c2ff6SSatish Balay        for (k=0; k<nzA; k++) {
1057233c2ff6SSatish Balay 	 ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1058233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1059233c2ff6SSatish Balay         }
1060233c2ff6SSatish Balay         for (k=0; k<nzB; k++) {
1061233c2ff6SSatish Balay 	  ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1062233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1063233c2ff6SSatish Balay         }
1064233c2ff6SSatish Balay #else
1065ca161407SBarry Smith         for (k=0; k<nzA; k++) {
106636f4e84dSSatish Balay           if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; }
1067ca161407SBarry Smith         }
1068ca161407SBarry Smith         for (k=0; k<nzB; k++) {
106936f4e84dSSatish Balay           if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; }
1070d5b485abSSatish Balay         }
1071233c2ff6SSatish Balay #endif
1072d5b485abSSatish Balay       }
1073a2feddc5SSatish Balay     }
1074ca161407SBarry Smith   }
1075233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1076d5b485abSSatish Balay   /* Create row map*/
1077b0a32e0cSBarry Smith   ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rowmaps);CHKERRQ(ierr);
1078ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
1079ff0794f7SSatish Balay     ierr = PetscTableCreate(nrow[i]+1,&rowmaps[i]);CHKERRQ(ierr);
1080233c2ff6SSatish Balay   }
1081233c2ff6SSatish Balay #else
1082233c2ff6SSatish Balay   /* Create row map*/
1083b24ad042SBarry Smith   len     = (1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt);
108482502324SSatish Balay   ierr    = PetscMalloc(len,&rmap);CHKERRQ(ierr);
1085b24ad042SBarry Smith   rmap[0] = (PetscInt*)(rmap + ismax);
1086b24ad042SBarry Smith   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr);
1087233c2ff6SSatish Balay   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;}
1088233c2ff6SSatish Balay #endif
1089d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1090d5b485abSSatish Balay     irow_i = irow[i];
1091d5b485abSSatish Balay     jmax   = nrow[i];
1092233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1093233c2ff6SSatish Balay     lrow1_grow1 = rowmaps[i];
1094233c2ff6SSatish Balay     for (j=0; j<jmax; j++) {
1095233c2ff6SSatish Balay       ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr);
1096233c2ff6SSatish Balay     }
1097233c2ff6SSatish Balay #else
1098233c2ff6SSatish Balay     rmap_i = rmap[i];
1099d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1100d5b485abSSatish Balay       rmap_i[irow_i[j]] = j;
1101d5b485abSSatish Balay     }
1102233c2ff6SSatish Balay #endif
1103d5b485abSSatish Balay   }
1104d5b485abSSatish Balay 
1105d5b485abSSatish Balay   /* Update lens from offproc data */
1106d5b485abSSatish Balay   {
1107b24ad042SBarry Smith     PetscInt    *rbuf2_i,*rbuf3_i,*sbuf1_i;
1108b24ad042SBarry Smith     PetscMPIInt ii;
1109d5b485abSSatish Balay 
1110d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1111b24ad042SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr);
1112b24ad042SBarry Smith       idex   = pa[ii];
1113999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1114d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1115d5b485abSSatish Balay       ct1     = 2*jmax+1;
1116d5b485abSSatish Balay       ct2     = 0;
1117b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1118b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
1119d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1120d5b485abSSatish Balay         is_no   = sbuf1_i[2*j-1];
1121d5b485abSSatish Balay         max1    = sbuf1_i[2*j];
1122d5b485abSSatish Balay         lens_i  = lens[is_no];
1123233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1124233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1125233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1126233c2ff6SSatish Balay #else
1127d5b485abSSatish Balay         cmap_i  = cmap[is_no];
1128d5b485abSSatish Balay 	rmap_i  = rmap[is_no];
1129233c2ff6SSatish Balay #endif
1130d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1131233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1132233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1133233c2ff6SSatish Balay           row--;
1134e005ede5SBarry Smith           if (row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); }
1135233c2ff6SSatish Balay #else
1136d5b485abSSatish Balay           row  = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
1137233c2ff6SSatish Balay #endif
1138d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1139d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1140233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1141233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
1142233c2ff6SSatish Balay 	    if (tt) {
1143233c2ff6SSatish Balay               lens_i[row]++;
1144233c2ff6SSatish Balay             }
1145233c2ff6SSatish Balay #else
1146d5b485abSSatish Balay             if (cmap_i[rbuf3_i[ct2]]) {
1147d5b485abSSatish Balay               lens_i[row]++;
1148d5b485abSSatish Balay             }
1149233c2ff6SSatish Balay #endif
1150d5b485abSSatish Balay           }
1151d5b485abSSatish Balay         }
1152d5b485abSSatish Balay       }
1153d5b485abSSatish Balay     }
1154d5b485abSSatish Balay   }
1155606d414cSSatish Balay   ierr = PetscFree(r_status3);CHKERRQ(ierr);
1156606d414cSSatish Balay   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
11570c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
1158606d414cSSatish Balay   ierr = PetscFree(s_status3);CHKERRQ(ierr);
1159606d414cSSatish Balay   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1160d5b485abSSatish Balay 
1161d5b485abSSatish Balay   /* Create the submatrices */
1162d5b485abSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1163d5b485abSSatish Balay     /*
1164d5b485abSSatish Balay         Assumes new rows are same length as the old rows, hence bug!
1165d5b485abSSatish Balay     */
1166d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1167df36731bSBarry Smith       mat = (Mat_SeqBAIJ *)(submats[i]->data);
1168d0f46423SBarry Smith       if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) {
116929bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1170d5b485abSSatish Balay       }
1171b24ad042SBarry Smith       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
1172abc0a331SBarry Smith       if (!flag) {
117329bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1174d5b485abSSatish Balay       }
1175d5b485abSSatish Balay       /* Initial matrix as if empty */
1176b24ad042SBarry Smith       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr);
1177ce60de66SLois Curfman McInnes       submats[i]->factor = C->factor;
1178d5b485abSSatish Balay     }
1179ca161407SBarry Smith   } else {
1180d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1181f69a0ea3SMatthew Knepley       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
1182f69a0ea3SMatthew Knepley       ierr = MatSetSizes(submats[i],nrow[i]*bs,ncol[i]*bs,nrow[i]*bs,ncol[i]*bs);CHKERRQ(ierr);
11837adad957SLisandro Dalcin       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
1184d0f46423SBarry Smith       ierr = MatSeqBAIJSetPreallocation(submats[i],C->rmap->bs,0,lens[i]);CHKERRQ(ierr);
1185d0f46423SBarry Smith       ierr = MatSeqSBAIJSetPreallocation(submats[i],C->rmap->bs,0,lens[i]);CHKERRQ(ierr);
1186d5b485abSSatish Balay     }
1187d5b485abSSatish Balay   }
1188d5b485abSSatish Balay 
1189d5b485abSSatish Balay   /* Assemble the matrices */
1190d5b485abSSatish Balay   /* First assemble the local rows */
1191d5b485abSSatish Balay   {
1192b24ad042SBarry Smith     PetscInt       ilen_row,*imat_ilen,*imat_j,*imat_i;
11933eda8832SBarry Smith     MatScalar *imat_a;
1194d5b485abSSatish Balay 
1195d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1196df36731bSBarry Smith       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1197d5b485abSSatish Balay       imat_ilen = mat->ilen;
1198d5b485abSSatish Balay       imat_j    = mat->j;
1199d5b485abSSatish Balay       imat_i    = mat->i;
1200d5b485abSSatish Balay       imat_a    = mat->a;
1201233c2ff6SSatish Balay 
1202233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1203233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
1204233c2ff6SSatish Balay       lrow1_grow1 = rowmaps[i];
1205233c2ff6SSatish Balay #else
1206d5b485abSSatish Balay       cmap_i    = cmap[i];
1207d5b485abSSatish Balay       rmap_i    = rmap[i];
1208233c2ff6SSatish Balay #endif
1209d5b485abSSatish Balay       irow_i    = irow[i];
1210d5b485abSSatish Balay       jmax      = nrow[i];
1211d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1212d5b485abSSatish Balay         row      = irow_i[j];
1213233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1214899cda47SBarry Smith 	ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1215233c2ff6SSatish Balay #else
1216d5b485abSSatish Balay 	proc = rtable[row];
1217233c2ff6SSatish Balay #endif
1218d5b485abSSatish Balay         if (proc == rank) {
121936f4e84dSSatish Balay           row      = row - rstart;
122036f4e84dSSatish Balay           nzA      = a_i[row+1] - a_i[row];
122136f4e84dSSatish Balay           nzB      = b_i[row+1] - b_i[row];
122236f4e84dSSatish Balay           cworkA   = a_j + a_i[row];
122336f4e84dSSatish Balay           cworkB   = b_j + b_i[row];
122436f4e84dSSatish Balay           vworkA   = a_a + a_i[row]*bs2;
122536f4e84dSSatish Balay           vworkB   = b_a + b_i[row]*bs2;
1226233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1227233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr);
1228233c2ff6SSatish Balay           row--;
1229e005ede5SBarry Smith           if (row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); }
1230233c2ff6SSatish Balay #else
123136f4e84dSSatish Balay           row      = rmap_i[row + rstart];
1232233c2ff6SSatish Balay #endif
1233df36731bSBarry Smith           mat_i    = imat_i[row];
123436f4e84dSSatish Balay           mat_a    = imat_a + mat_i*bs2;
1235d5b485abSSatish Balay           mat_j    = imat_j + mat_i;
123636f4e84dSSatish Balay           ilen_row = imat_ilen[row];
123736f4e84dSSatish Balay 
123836f4e84dSSatish Balay           /* load the column indices for this row into cols*/
123936f4e84dSSatish Balay           for (l=0; l<nzB; l++) {
124036f4e84dSSatish Balay 	    if ((ctmp = bmap[cworkB[l]]) < cstart) {
1241233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1242233c2ff6SSatish Balay 	      ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr);
1243233c2ff6SSatish Balay 	      if (tcol) {
1244233c2ff6SSatish Balay #else
124536f4e84dSSatish Balay               if ((tcol = cmap_i[ctmp])) {
1246233c2ff6SSatish Balay #endif
1247df36731bSBarry Smith                 *mat_j++ = tcol - 1;
12483eda8832SBarry Smith                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1249549d3d68SSatish Balay                 mat_a   += bs2;
1250d5b485abSSatish Balay                 ilen_row++;
1251d5b485abSSatish Balay               }
1252ca161407SBarry Smith             } else break;
125336f4e84dSSatish Balay           }
125436f4e84dSSatish Balay           imark = l;
125536f4e84dSSatish Balay           for (l=0; l<nzA; l++) {
1256233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1257233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1258233c2ff6SSatish Balay 	    if (tcol) {
1259233c2ff6SSatish Balay #else
126036f4e84dSSatish Balay 	    if ((tcol = cmap_i[cstart + cworkA[l]])) {
1261233c2ff6SSatish Balay #endif
126236f4e84dSSatish Balay               *mat_j++ = tcol - 1;
12633eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1264549d3d68SSatish Balay               mat_a   += bs2;
126536f4e84dSSatish Balay               ilen_row++;
126636f4e84dSSatish Balay             }
126736f4e84dSSatish Balay           }
126836f4e84dSSatish Balay           for (l=imark; l<nzB; l++) {
1269233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1270233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1271233c2ff6SSatish Balay 	    if (tcol) {
1272233c2ff6SSatish Balay #else
127336f4e84dSSatish Balay             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1274233c2ff6SSatish Balay #endif
127536f4e84dSSatish Balay               *mat_j++ = tcol - 1;
12763eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1277549d3d68SSatish Balay               mat_a   += bs2;
127836f4e84dSSatish Balay               ilen_row++;
127936f4e84dSSatish Balay             }
128036f4e84dSSatish Balay           }
1281d5b485abSSatish Balay           imat_ilen[row] = ilen_row;
1282d5b485abSSatish Balay         }
1283d5b485abSSatish Balay       }
128436f4e84dSSatish Balay 
1285d5b485abSSatish Balay     }
1286d5b485abSSatish Balay   }
1287d5b485abSSatish Balay 
1288d5b485abSSatish Balay   /*   Now assemble the off proc rows*/
1289d5b485abSSatish Balay   {
1290b24ad042SBarry Smith     PetscInt    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1291b24ad042SBarry Smith     PetscInt    *imat_j,*imat_i;
12923eda8832SBarry Smith     MatScalar   *imat_a,*rbuf4_i;
1293b24ad042SBarry Smith     PetscMPIInt ii;
1294d5b485abSSatish Balay 
1295d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1296b24ad042SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr);
1297b24ad042SBarry Smith       idex   = pa[ii];
1298999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1299d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1300d5b485abSSatish Balay       ct1     = 2*jmax + 1;
1301d5b485abSSatish Balay       ct2     = 0;
1302b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1303b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
1304b24ad042SBarry Smith       rbuf4_i = rbuf4[ii];
1305d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1306d5b485abSSatish Balay         is_no     = sbuf1_i[2*j-1];
1307233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1308233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1309233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1310233c2ff6SSatish Balay #else
1311d5b485abSSatish Balay         rmap_i    = rmap[is_no];
1312d5b485abSSatish Balay         cmap_i    = cmap[is_no];
1313233c2ff6SSatish Balay #endif
1314df36731bSBarry Smith         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1315d5b485abSSatish Balay         imat_ilen = mat->ilen;
1316d5b485abSSatish Balay         imat_j    = mat->j;
1317d5b485abSSatish Balay         imat_i    = mat->i;
1318d5b485abSSatish Balay         imat_a    = mat->a;
1319d5b485abSSatish Balay         max1      = sbuf1_i[2*j];
1320d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1321d5b485abSSatish Balay           row   = sbuf1_i[ct1];
1322233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1323233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr);
1324233c2ff6SSatish Balay           row--;
1325e005ede5SBarry Smith           if(row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); }
1326233c2ff6SSatish Balay #else
1327d5b485abSSatish Balay           row   = rmap_i[row];
1328233c2ff6SSatish Balay #endif
1329d5b485abSSatish Balay           ilen  = imat_ilen[row];
1330df36731bSBarry Smith           mat_i = imat_i[row];
133136f4e84dSSatish Balay           mat_a = imat_a + mat_i*bs2;
1332d5b485abSSatish Balay           mat_j = imat_j + mat_i;
1333d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1334d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1335233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1336233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1337233c2ff6SSatish Balay 	    if (tcol) {
1338233c2ff6SSatish Balay #else
1339d5b485abSSatish Balay 	    if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1340233c2ff6SSatish Balay #endif
1341df36731bSBarry Smith               *mat_j++    = tcol - 1;
134236f4e84dSSatish Balay               /* *mat_a++= rbuf4_i[ct2]; */
13433eda8832SBarry Smith               ierr        = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1344549d3d68SSatish Balay               mat_a      += bs2;
1345d5b485abSSatish Balay               ilen++;
1346d5b485abSSatish Balay             }
1347d5b485abSSatish Balay           }
1348d5b485abSSatish Balay           imat_ilen[row] = ilen;
1349d5b485abSSatish Balay         }
1350d5b485abSSatish Balay       }
1351d5b485abSSatish Balay     }
1352d5b485abSSatish Balay   }
1353606d414cSSatish Balay   ierr = PetscFree(r_status4);CHKERRQ(ierr);
1354606d414cSSatish Balay   ierr = PetscFree(r_waits4);CHKERRQ(ierr);
13550c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
1356606d414cSSatish Balay   ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1357606d414cSSatish Balay   ierr = PetscFree(s_status4);CHKERRQ(ierr);
1358d5b485abSSatish Balay 
1359d5b485abSSatish Balay   /* Restore the indices */
1360d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1361d5b485abSSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1362d5b485abSSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1363d5b485abSSatish Balay   }
1364d5b485abSSatish Balay 
1365d5b485abSSatish Balay   /* Destroy allocated memory */
1366606d414cSSatish Balay   ierr = PetscFree(irow);CHKERRQ(ierr);
1367606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
1368606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
1369d5b485abSSatish Balay 
1370606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1371606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1372d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
1373606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1374d5b485abSSatish Balay   }
1375d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1376606d414cSSatish Balay     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1377606d414cSSatish Balay     ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);
1378d5b485abSSatish Balay   }
1379d5b485abSSatish Balay 
1380606d414cSSatish Balay   ierr = PetscFree(sbuf2);CHKERRQ(ierr);
1381606d414cSSatish Balay   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1382606d414cSSatish Balay   ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1383606d414cSSatish Balay   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1384606d414cSSatish Balay   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
1385606d414cSSatish Balay   ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1386606d414cSSatish Balay   ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
1387d5b485abSSatish Balay 
1388233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1389ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
13909c666560SBarry Smith     ierr = PetscTableDestroy(rowmaps[i]);CHKERRQ(ierr);
13919c666560SBarry Smith     ierr = PetscTableDestroy(colmaps[i]);CHKERRQ(ierr);
1392233c2ff6SSatish Balay   }
1393233c2ff6SSatish Balay   ierr = PetscFree(colmaps);CHKERRQ(ierr);
1394233c2ff6SSatish Balay   ierr = PetscFree(rowmaps);CHKERRQ(ierr);
1395233c2ff6SSatish Balay #else
1396606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
1397233c2ff6SSatish Balay   ierr = PetscFree(cmap);CHKERRQ(ierr);
1398233c2ff6SSatish Balay #endif
1399606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1400d5b485abSSatish Balay 
1401d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
140236f4e84dSSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
140336f4e84dSSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1404d5b485abSSatish Balay   }
14053a40ed3dSBarry Smith   PetscFunctionReturn(0);
1406d5b485abSSatish Balay }
1407ca161407SBarry Smith 
1408