xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision c05d87d6c2fa1d0952fe676de20339d09daa79bc)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
323bdbc58SBarry Smith 
4d5b485abSSatish Balay /*
5d5b485abSSatish Balay    Routines to compute overlapping regions of a parallel MPI matrix
6d5b485abSSatish Balay   and to find submatrices that were shared across processors.
7d5b485abSSatish Balay */
87c4f633dSBarry Smith #include "../src/mat/impls/baij/mpi/mpibaij.h"
90a835dfdSSatish Balay #include "petscbt.h"
10d5b485abSSatish Balay 
11b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat,PetscInt,IS *);
12b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char **,PetscInt*,PetscInt**);
13b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt **,PetscInt**,PetscInt*);
14b24ad042SBarry Smith EXTERN PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
15b24ad042SBarry Smith EXTERN PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
16d5b485abSSatish Balay 
174a2ae208SSatish Balay #undef __FUNCT__
184a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ"
19b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov)
20d5b485abSSatish Balay {
216849ba73SBarry Smith   PetscErrorCode ierr;
22d0f46423SBarry Smith   PetscInt       i,N=C->cmap->N, bs=C->rmap->bs;
2336f4e84dSSatish Balay   IS             *is_new;
2436f4e84dSSatish Balay 
253a40ed3dSBarry Smith   PetscFunctionBegin;
2682502324SSatish Balay   ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr);
27df36731bSBarry Smith   /* Convert the indices into block format */
2827f478b1SHong Zhang   ierr = ISCompressIndicesGeneral(N,bs,imax,is,is_new);CHKERRQ(ierr);
2929bbc08cSBarry Smith   if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
30d5b485abSSatish Balay   for (i=0; i<ov; ++i) {
3136f4e84dSSatish Balay     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr);
32d5b485abSSatish Balay   }
33ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);}
3427f478b1SHong Zhang   ierr = ISExpandIndicesGeneral(N,bs,imax,is_new,is);CHKERRQ(ierr);
35ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);}
36606d414cSSatish Balay   ierr = PetscFree(is_new);CHKERRQ(ierr);
373a40ed3dSBarry Smith   PetscFunctionReturn(0);
38d5b485abSSatish Balay }
39d5b485abSSatish Balay 
40d5b485abSSatish Balay /*
41d5b485abSSatish Balay   Sample message format:
42d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
430e9b0e7eSHong Zhang   to index sets is[1], is[5]
44d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
45d5b485abSSatish Balay   -----------
46d5b485abSSatish Balay   mesg [1] = 1 => is[1]
47d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
48d5b485abSSatish Balay   -----------
49d5b485abSSatish Balay   mesg [5] = 5  => is[5]
50d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
51d5b485abSSatish Balay   -----------
52d5b485abSSatish Balay   mesg [7]
530e9b0e7eSHong Zhang   mesg [n]  data(is[1])
54d5b485abSSatish Balay   -----------
55d5b485abSSatish Balay   mesg[n+1]
56d5b485abSSatish Balay   mesg[m]  data(is[5])
57d5b485abSSatish Balay   -----------
58d5b485abSSatish Balay 
59d5b485abSSatish Balay   Notes:
60d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
61d5b485abSSatish Balay   nrqr - no of requests recieved (which have to be or which have been processed
62d5b485abSSatish Balay */
634a2ae208SSatish Balay #undef __FUNCT__
644a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once"
65b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[])
66d5b485abSSatish Balay {
67df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
685d0c19d7SBarry Smith   const PetscInt **idx,*idx_i;
695d0c19d7SBarry Smith   PetscInt       *n,*w3,*w4,*rtable,**data,len;
706849ba73SBarry Smith   PetscErrorCode ierr;
71b24ad042SBarry Smith   PetscMPIInt    size,rank,tag1,tag2,*w2,*w1,nrqr;
72b24ad042SBarry Smith   PetscInt       Mbs,i,j,k,**rbuf,row,proc,nrqs,msz,**outdat,**ptr;
73b24ad042SBarry Smith   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2;
74b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1,*onodes2,*olengths2;
75f1af5d2fSBarry Smith   PetscBT        *table;
76d5b485abSSatish Balay   MPI_Comm       comm;
77d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
78d5b485abSSatish Balay   MPI_Status     *s_status,*recv_status;
79d5b485abSSatish Balay 
803a40ed3dSBarry Smith   PetscFunctionBegin;
817adad957SLisandro Dalcin   comm   = ((PetscObject)C)->comm;
82d5b485abSSatish Balay   size   = c->size;
83d5b485abSSatish Balay   rank   = c->rank;
84df36731bSBarry Smith   Mbs    = c->Mbs;
85d5b485abSSatish Balay 
86c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
87c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
88c7dd2924SSatish Balay 
8923bdbc58SBarry Smith   ierr = PetscMalloc3(imax+1,const PetscInt*,&idx,imax,PetscInt,&n,Mbs,PetscInt,&rtable);CHKERRQ(ierr);
90d5b485abSSatish Balay   for (i=0; i<imax; i++) {
91d5b485abSSatish Balay     ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr);
92b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
93d5b485abSSatish Balay   }
94d5b485abSSatish Balay 
95d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
96d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
97899cda47SBarry Smith     len = c->rangebs[i+1];
98d5b485abSSatish Balay     for (; j<len; j++) {
99d5b485abSSatish Balay       rtable[j] = i;
100d5b485abSSatish Balay     }
101d5b485abSSatish Balay   }
102d5b485abSSatish Balay 
103d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
104d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
10523bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr);
10623bdbc58SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
10723bdbc58SBarry Smith   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
10823bdbc58SBarry Smith   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
109d5b485abSSatish Balay   for (i=0; i<imax; i++) {
110b24ad042SBarry Smith     ierr  = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
111d5b485abSSatish Balay     idx_i = idx[i];
112d5b485abSSatish Balay     len   = n[i];
113d5b485abSSatish Balay     for (j=0; j<len; j++) {
114d5b485abSSatish Balay       row  = idx_i[j];
1156b41c64dSBarry Smith       if (row < 0) {
116e005ede5SBarry Smith         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries");
1176b41c64dSBarry Smith       }
118d5b485abSSatish Balay       proc = rtable[row];
119d5b485abSSatish Balay       w4[proc]++;
120d5b485abSSatish Balay     }
121d5b485abSSatish Balay     for (j=0; j<size; j++){
122d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
123d5b485abSSatish Balay     }
124d5b485abSSatish Balay   }
125d5b485abSSatish Balay 
126d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
127d5b485abSSatish Balay   msz      = 0;              /* total mesg length (for all proc */
1280e9b0e7eSHong Zhang   w1[rank] = 0;              /* no mesg sent to itself */
129d5b485abSSatish Balay   w3[rank] = 0;
130d5b485abSSatish Balay   for (i=0; i<size; i++) {
131d5b485abSSatish Balay     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
132d5b485abSSatish Balay   }
133d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
134b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr);
135d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
136d5b485abSSatish Balay     if (w1[i]) {pa[j] = i; j++;}
137d5b485abSSatish Balay   }
138d5b485abSSatish Balay 
139d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
140d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
141d5b485abSSatish Balay     j      = pa[i];
142d5b485abSSatish Balay     w1[j] += w2[j] + 2*w3[j];
143d5b485abSSatish Balay     msz   += w1[j];
144d5b485abSSatish Balay   }
145d5b485abSSatish Balay 
146c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
147a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
148c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
149d5b485abSSatish Balay 
150c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
151c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr);
152d5b485abSSatish Balay 
153d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
15423bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscInt*,&outdat,size,PetscInt*,&ptr,msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr);
15523bdbc58SBarry Smith   ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr);
15623bdbc58SBarry Smith   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
157d5b485abSSatish Balay   {
158b24ad042SBarry Smith     PetscInt *iptr = tmp,ict  = 0;
159d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
160d5b485abSSatish Balay       j         = pa[i];
161d5b485abSSatish Balay       iptr     +=  ict;
162d5b485abSSatish Balay       outdat[j] = iptr;
163d5b485abSSatish Balay       ict       = w1[j];
164d5b485abSSatish Balay     }
165d5b485abSSatish Balay   }
166d5b485abSSatish Balay 
167d5b485abSSatish Balay   /* Form the outgoing messages */
168d5b485abSSatish Balay   /*plug in the headers*/
169d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
170d5b485abSSatish Balay     j            = pa[i];
171d5b485abSSatish Balay     outdat[j][0] = 0;
172b24ad042SBarry Smith     ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
173d5b485abSSatish Balay     ptr[j]       = outdat[j] + 2*w3[j] + 1;
174d5b485abSSatish Balay   }
175d5b485abSSatish Balay 
176d5b485abSSatish Balay   /* Memory for doing local proc's work*/
177d5b485abSSatish Balay   {
178b24ad042SBarry Smith     PetscInt  *d_p;
179d5b485abSSatish Balay     char      *t_p;
180d5b485abSSatish Balay 
18123bdbc58SBarry Smith     /* should change to use PetscMallocN() */
182052f0c41SBarry Smith     ierr = PetscMalloc((imax)*(sizeof(PetscBT) + sizeof(PetscInt*)+ sizeof(PetscInt)) +
183052f0c41SBarry Smith       (Mbs)*imax*sizeof(PetscInt)  + (Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char),&table);CHKERRQ(ierr);
184052f0c41SBarry Smith     ierr = PetscMemzero(table,(imax)*(sizeof(PetscBT) + sizeof(PetscInt*)+ sizeof(PetscInt)) +
185052f0c41SBarry Smith       (Mbs)*imax*sizeof(PetscInt)  + (Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char));CHKERRQ(ierr);
186b24ad042SBarry Smith     data = (PetscInt **)(table + imax);
187b24ad042SBarry Smith     isz  = (PetscInt  *)(data  + imax);
188b24ad042SBarry Smith     d_p  = (PetscInt  *)(isz   + imax);
189bbd702dbSSatish Balay     t_p  = (char *)(d_p   + Mbs*imax);
190bbd702dbSSatish Balay     for (i=0; i<imax; i++) {
191b6410449SSatish Balay       table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
192bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
193d5b485abSSatish Balay     }
194d5b485abSSatish Balay   }
195d5b485abSSatish Balay 
196d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
197d5b485abSSatish Balay   {
198b24ad042SBarry Smith     PetscInt     n_i,*data_i,isz_i,*outdat_j,ctr_j;
199f1af5d2fSBarry Smith     PetscBT table_i;
200d5b485abSSatish Balay 
201d5b485abSSatish Balay     for (i=0; i<imax; i++) {
202b24ad042SBarry Smith       ierr    = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
203d5b485abSSatish Balay       n_i     = n[i];
204d5b485abSSatish Balay       table_i = table[i];
205d5b485abSSatish Balay       idx_i   = idx[i];
206d5b485abSSatish Balay       data_i  = data[i];
207d5b485abSSatish Balay       isz_i   = isz[i];
208d5b485abSSatish Balay       for (j=0;  j<n_i; j++) {  /* parse the indices of each IS */
209d5b485abSSatish Balay         row  = idx_i[j];
210d5b485abSSatish Balay         proc = rtable[row];
211d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
212d5b485abSSatish Balay           ctr[proc]++;
213d5b485abSSatish Balay           *ptr[proc] = row;
214d5b485abSSatish Balay           ptr[proc]++;
215d5b485abSSatish Balay         }
216d5b485abSSatish Balay         else { /* Update the local table */
217f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
218d5b485abSSatish Balay         }
219d5b485abSSatish Balay       }
220d5b485abSSatish Balay       /* Update the headers for the current IS */
221d5b485abSSatish Balay       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
222d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
223d5b485abSSatish Balay           outdat_j        = outdat[j];
224d5b485abSSatish Balay           k               = ++outdat_j[0];
225d5b485abSSatish Balay           outdat_j[2*k]   = ctr_j;
226d5b485abSSatish Balay           outdat_j[2*k-1] = i;
227d5b485abSSatish Balay         }
228d5b485abSSatish Balay       }
229d5b485abSSatish Balay       isz[i] = isz_i;
230d5b485abSSatish Balay     }
231d5b485abSSatish Balay   }
232d5b485abSSatish Balay 
233d5b485abSSatish Balay   /*  Now  post the sends */
234b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
235d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
236d5b485abSSatish Balay     j    = pa[i];
237b24ad042SBarry Smith     ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr);
238d5b485abSSatish Balay   }
239d5b485abSSatish Balay 
240d5b485abSSatish Balay   /* No longer need the original indices*/
241d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
242d5b485abSSatish Balay     ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr);
243d5b485abSSatish Balay   }
24423bdbc58SBarry Smith   ierr = PetscFree3(idx,n,rtable);CHKERRQ(ierr);
245d5b485abSSatish Balay 
246d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
247d5b485abSSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
248d5b485abSSatish Balay   }
249d5b485abSSatish Balay 
250d5b485abSSatish Balay   /* Do Local work*/
251df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
252d5b485abSSatish Balay 
253d5b485abSSatish Balay   /* Receive messages*/
254b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr);
2550c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);}
256d5b485abSSatish Balay 
257b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr);
2580c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);}
259d5b485abSSatish Balay 
260d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
26123bdbc58SBarry Smith   ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr);
26223bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
263d5b485abSSatish Balay 
264b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&xdata);CHKERRQ(ierr);
265b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt),&isz1);CHKERRQ(ierr);
266df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
267*c05d87d6SBarry Smith   ierr = PetscFree(rbuf[0]);CHKERRQ(ierr);
268606d414cSSatish Balay   ierr = PetscFree(rbuf);CHKERRQ(ierr);
269d5b485abSSatish Balay 
270d5b485abSSatish Balay   /* Send the data back*/
271d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
272d5b485abSSatish Balay   {
273b24ad042SBarry Smith     PetscMPIInt *rw1;
274d5b485abSSatish Balay 
275b24ad042SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&rw1);CHKERRQ(ierr);
276b24ad042SBarry Smith     ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr);
277c7dd2924SSatish Balay 
278d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
279d5b485abSSatish Balay       proc      = recv_status[i].MPI_SOURCE;
280e005ede5SBarry Smith       if (proc != onodes1[i]) SETERRQ(PETSC_ERR_PLIB,"MPI_SOURCE mismatch");
281d5b485abSSatish Balay       rw1[proc] = isz1[i];
282d5b485abSSatish Balay     }
283d5b485abSSatish Balay 
284c7dd2924SSatish Balay     ierr = PetscFree(onodes1);CHKERRQ(ierr);
285c7dd2924SSatish Balay     ierr = PetscFree(olengths1);CHKERRQ(ierr);
286d5b485abSSatish Balay 
287c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
288c7dd2924SSatish Balay     ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr);
28903512dc5SSatish Balay     ierr = PetscFree(rw1);CHKERRQ(ierr);
290c7dd2924SSatish Balay   }
291c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
292c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr);
293d5b485abSSatish Balay 
294d5b485abSSatish Balay   /*  Now  post the sends */
295b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
296d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
297d5b485abSSatish Balay     j    = recv_status[i].MPI_SOURCE;
298b24ad042SBarry Smith     ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr);
299d5b485abSSatish Balay   }
300d5b485abSSatish Balay 
301d5b485abSSatish Balay   /* receive work done on other processors*/
302d5b485abSSatish Balay   {
303b24ad042SBarry Smith     PetscMPIInt idex;
304b24ad042SBarry Smith     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
305f1af5d2fSBarry Smith     PetscBT     table_i;
306d5b485abSSatish Balay     MPI_Status  *status2;
307d5b485abSSatish Balay 
308169449a0SSatish Balay     ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr);
309d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
310999d9058SBarry Smith       ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr);
311d5b485abSSatish Balay       /* Process the message*/
312999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
313d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
314999d9058SBarry Smith       jmax    = rbuf2[idex][0];
315d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
316d5b485abSSatish Balay         max     = rbuf2_i[2*j];
317d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
318d5b485abSSatish Balay         isz_i   = isz[is_no];
319d5b485abSSatish Balay         data_i  = data[is_no];
320d5b485abSSatish Balay         table_i = table[is_no];
321d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
322d5b485abSSatish Balay           row = rbuf2_i[ct1];
323f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
324d5b485abSSatish Balay         }
325d5b485abSSatish Balay         isz[is_no] = isz_i;
326d5b485abSSatish Balay       }
327d5b485abSSatish Balay     }
3280c468ba9SBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);}
329606d414cSSatish Balay     ierr = PetscFree(status2);CHKERRQ(ierr);
330d5b485abSSatish Balay   }
331d5b485abSSatish Balay 
332d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
333029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr);
334d5b485abSSatish Balay   }
335d5b485abSSatish Balay 
336c7dd2924SSatish Balay 
337c7dd2924SSatish Balay   ierr = PetscFree(onodes2);CHKERRQ(ierr);
338c7dd2924SSatish Balay   ierr = PetscFree(olengths2);CHKERRQ(ierr);
339c7dd2924SSatish Balay 
340606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
341*c05d87d6SBarry Smith   ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr);
342606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
343606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
344606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
345606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
346606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
347606d414cSSatish Balay   ierr = PetscFree(table);CHKERRQ(ierr);
348606d414cSSatish Balay   ierr = PetscFree(s_status);CHKERRQ(ierr);
349606d414cSSatish Balay   ierr = PetscFree(recv_status);CHKERRQ(ierr);
350606d414cSSatish Balay   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
351606d414cSSatish Balay   ierr = PetscFree(xdata);CHKERRQ(ierr);
352606d414cSSatish Balay   ierr = PetscFree(isz1);CHKERRQ(ierr);
3533a40ed3dSBarry Smith   PetscFunctionReturn(0);
354d5b485abSSatish Balay }
355d5b485abSSatish Balay 
3564a2ae208SSatish Balay #undef __FUNCT__
3574a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local"
358d5b485abSSatish Balay /*
359df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
360d5b485abSSatish Balay        the work on the local processor.
361d5b485abSSatish Balay 
362d5b485abSSatish Balay      Inputs:
363df36731bSBarry Smith       C      - MAT_MPIBAIJ;
364d5b485abSSatish Balay       imax - total no of index sets processed at a time;
365df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
366d5b485abSSatish Balay 
367d5b485abSSatish Balay      Output:
3680e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
369d5b485abSSatish Balay                to each index set;
370d5b485abSSatish Balay       data   - pointer to the solutions
371d5b485abSSatish Balay */
372b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data)
373d5b485abSSatish Balay {
374df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
375d5b485abSSatish Balay   Mat         A = c->A,B = c->B;
376df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
377b24ad042SBarry Smith   PetscInt    start,end,val,max,rstart,cstart,*ai,*aj;
378b24ad042SBarry Smith   PetscInt    *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
379f1af5d2fSBarry Smith   PetscBT     table_i;
380d5b485abSSatish Balay 
3813a40ed3dSBarry Smith   PetscFunctionBegin;
382899cda47SBarry Smith   rstart = c->rstartbs;
383899cda47SBarry Smith   cstart = c->cstartbs;
384d5b485abSSatish Balay   ai     = a->i;
385df36731bSBarry Smith   aj     = a->j;
386d5b485abSSatish Balay   bi     = b->i;
387df36731bSBarry Smith   bj     = b->j;
388d5b485abSSatish Balay   garray = c->garray;
389d5b485abSSatish Balay 
390d5b485abSSatish Balay 
391d5b485abSSatish Balay   for (i=0; i<imax; i++) {
392d5b485abSSatish Balay     data_i  = data[i];
393d5b485abSSatish Balay     table_i = table[i];
394d5b485abSSatish Balay     isz_i   = isz[i];
395d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
396d5b485abSSatish Balay       row   = data_i[j] - rstart;
397d5b485abSSatish Balay       start = ai[row];
398d5b485abSSatish Balay       end   = ai[row+1];
399d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
400df36731bSBarry Smith         val = aj[k] + cstart;
401f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
402d5b485abSSatish Balay       }
403d5b485abSSatish Balay       start = bi[row];
404d5b485abSSatish Balay       end   = bi[row+1];
405d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
406df36731bSBarry Smith         val = garray[bj[k]];
407f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
408d5b485abSSatish Balay       }
409d5b485abSSatish Balay     }
410d5b485abSSatish Balay     isz[i] = isz_i;
411d5b485abSSatish Balay   }
4123a40ed3dSBarry Smith   PetscFunctionReturn(0);
413d5b485abSSatish Balay }
4144a2ae208SSatish Balay #undef __FUNCT__
4154a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive"
416d5b485abSSatish Balay /*
417df36731bSBarry Smith       MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages,
418d5b485abSSatish Balay          and return the output
419d5b485abSSatish Balay 
420d5b485abSSatish Balay          Input:
421d5b485abSSatish Balay            C    - the matrix
422d5b485abSSatish Balay            nrqr - no of messages being processed.
423d5b485abSSatish Balay            rbuf - an array of pointers to the recieved requests
424d5b485abSSatish Balay 
425d5b485abSSatish Balay          Output:
426d5b485abSSatish Balay            xdata - array of messages to be sent back
427d5b485abSSatish Balay            isz1  - size of each message
428d5b485abSSatish Balay 
429a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
430d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
4310e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of
432d5b485abSSatish Balay memory is used.
433d5b485abSSatish Balay 
434d5b485abSSatish Balay */
435b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
436d5b485abSSatish Balay {
437df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
438d5b485abSSatish Balay   Mat            A = c->A,B = c->B;
439df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
4406849ba73SBarry Smith   PetscErrorCode ierr;
441b24ad042SBarry Smith   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
442b24ad042SBarry Smith   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
443b24ad042SBarry Smith   PetscInt       val,max1,max2,rank,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
444b24ad042SBarry Smith   PetscInt       *rbuf_i,kmax,rbuf_0;
445f1af5d2fSBarry Smith   PetscBT        xtable;
446d5b485abSSatish Balay 
4473a40ed3dSBarry Smith   PetscFunctionBegin;
448d5b485abSSatish Balay   rank   = c->rank;
449df36731bSBarry Smith   Mbs    = c->Mbs;
450899cda47SBarry Smith   rstart = c->rstartbs;
451899cda47SBarry Smith   cstart = c->cstartbs;
452d5b485abSSatish Balay   ai     = a->i;
453df36731bSBarry Smith   aj     = a->j;
454d5b485abSSatish Balay   bi     = b->i;
455df36731bSBarry Smith   bj     = b->j;
456d5b485abSSatish Balay   garray = c->garray;
457d5b485abSSatish Balay 
458d5b485abSSatish Balay 
459d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
460d5b485abSSatish Balay     rbuf_i  =  rbuf[i];
461d5b485abSSatish Balay     rbuf_0  =  rbuf_i[0];
462d5b485abSSatish Balay     ct     += rbuf_0;
463d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; }
464d5b485abSSatish Balay   }
465d5b485abSSatish Balay 
466701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
467701b8814SBarry Smith   else        max1 = 1;
468d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
469b24ad042SBarry Smith   ierr         = PetscMalloc(mem_estimate*sizeof(PetscInt),&xdata[0]);CHKERRQ(ierr);
470d5b485abSSatish Balay   ++no_malloc;
4716831982aSBarry Smith   ierr         = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr);
472b24ad042SBarry Smith   ierr         = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr);
473d5b485abSSatish Balay 
474d5b485abSSatish Balay   ct3 = 0;
475d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
476d5b485abSSatish Balay     rbuf_i =  rbuf[i];
477d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
478d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
479d5b485abSSatish Balay     ct2    =  ct1;
480d5b485abSSatish Balay     ct3    += ct1;
481d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
4826831982aSBarry Smith       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
483d5b485abSSatish Balay       oct2 = ct2;
484d5b485abSSatish Balay       kmax = rbuf_i[2*j];
485d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
486d5b485abSSatish Balay         row = rbuf_i[ct1];
487f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
488d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
489b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
490b24ad042SBarry Smith             ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
491b24ad042SBarry Smith             ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
492606d414cSSatish Balay             ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
493d5b485abSSatish Balay             xdata[0]     = tmp;
494d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
495d5b485abSSatish Balay             for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
496d5b485abSSatish Balay           }
497d5b485abSSatish Balay           xdata[i][ct2++] = row;
498d5b485abSSatish Balay           ct3++;
499d5b485abSSatish Balay         }
500d5b485abSSatish Balay       }
501d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
502d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
503d5b485abSSatish Balay         start = ai[row];
504d5b485abSSatish Balay         end   = ai[row+1];
505d5b485abSSatish Balay         for (l=start; l<end; l++) {
506df36731bSBarry Smith           val = aj[l] + cstart;
507f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
508d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
509b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
510b24ad042SBarry Smith               ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
511b24ad042SBarry Smith               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
512606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
513d5b485abSSatish Balay               xdata[0]     = tmp;
514d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
515d5b485abSSatish Balay               for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
516d5b485abSSatish Balay             }
517d5b485abSSatish Balay             xdata[i][ct2++] = val;
518d5b485abSSatish Balay             ct3++;
519d5b485abSSatish Balay           }
520d5b485abSSatish Balay         }
521d5b485abSSatish Balay         start = bi[row];
522d5b485abSSatish Balay         end   = bi[row+1];
523d5b485abSSatish Balay         for (l=start; l<end; l++) {
524df36731bSBarry Smith           val = garray[bj[l]];
525f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
526d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
527b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
528b24ad042SBarry Smith               ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
529b24ad042SBarry Smith               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
530606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
531d5b485abSSatish Balay               xdata[0]     = tmp;
532d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
533d5b485abSSatish Balay               for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
534d5b485abSSatish Balay             }
535d5b485abSSatish Balay             xdata[i][ct2++] = val;
536d5b485abSSatish Balay             ct3++;
537d5b485abSSatish Balay           }
538d5b485abSSatish Balay         }
539d5b485abSSatish Balay       }
540d5b485abSSatish Balay       /* Update the header*/
541d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
542d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
543d5b485abSSatish Balay     }
544d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
545d5b485abSSatish Balay     xdata[i+1]  = xdata[i] + ct2;
546d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
547d5b485abSSatish Balay   }
5486831982aSBarry Smith   ierr = PetscBTDestroy(xtable);CHKERRQ(ierr);
5491e2582c4SBarry Smith   ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr);
5503a40ed3dSBarry Smith   PetscFunctionReturn(0);
551d5b485abSSatish Balay }
552d5b485abSSatish Balay 
553b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *);
554a2feddc5SSatish Balay 
5554a2ae208SSatish Balay #undef __FUNCT__
5564a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ"
557b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
558d5b485abSSatish Balay {
55936f4e84dSSatish Balay   IS             *isrow_new,*iscol_new;
560cf4f527aSSatish Balay   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
5616849ba73SBarry Smith   PetscErrorCode ierr;
562d0f46423SBarry Smith   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,N=C->cmap->N,bs=C->rmap->bs;
563a2feddc5SSatish Balay 
5643a40ed3dSBarry Smith   PetscFunctionBegin;
5653a40ed3dSBarry Smith   /* The compression and expansion should be avoided. Does'nt point
566a2feddc5SSatish Balay      out errors might change the indices hence buggey */
567a2feddc5SSatish Balay 
56823bdbc58SBarry Smith   ierr = PetscMalloc2(ismax+1,IS,&isrow_new,ismax+1,IS,&iscol_new);CHKERRQ(ierr);
569d4503c43SBarry Smith   ierr = ISCompressIndicesGeneral(N,bs,ismax,isrow,isrow_new);CHKERRQ(ierr);
570d4503c43SBarry Smith   ierr = ISCompressIndicesGeneral(N,bs,ismax,iscol,iscol_new);CHKERRQ(ierr);
571cf4f527aSSatish Balay 
572cf4f527aSSatish Balay   /* Allocate memory to hold all the submatrices */
573cf4f527aSSatish Balay   if (scall != MAT_REUSE_MATRIX) {
57482502324SSatish Balay     ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr);
575cf4f527aSSatish Balay   }
576cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
577b24ad042SBarry Smith   nmax          = 20*1000000 / (c->Nbs * sizeof(PetscInt));
578329f5518SBarry Smith   if (!nmax) nmax = 1;
579cf4f527aSSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
580cf4f527aSSatish Balay 
581653e4784SBarry Smith   /* Make sure every processor loops through the nstages */
5827adad957SLisandro Dalcin   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);CHKERRQ(ierr);
583cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
584cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
585cf4f527aSSatish Balay     else if (pos == ismax) max_no = 0;
586cf4f527aSSatish Balay     else                   max_no = ismax-pos;
587cf4f527aSSatish Balay     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr);
588cf4f527aSSatish Balay     pos += max_no;
589cf4f527aSSatish Balay   }
59036f4e84dSSatish Balay 
59136f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
592ca161407SBarry Smith     ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr);
593ca161407SBarry Smith     ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr);
59436f4e84dSSatish Balay   }
59523bdbc58SBarry Smith   ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr);
5963a40ed3dSBarry Smith   PetscFunctionReturn(0);
597a2feddc5SSatish Balay }
598a2feddc5SSatish Balay 
599233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
6004a2ae208SSatish Balay #undef __FUNCT__
6014a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc"
602e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
603233c2ff6SSatish Balay {
604e005ede5SBarry Smith   PetscInt    nGlobalNd = proc_gnode[size];
605b9f7ace7SBarry Smith   PetscMPIInt fproc = PetscMPIIntCast( (PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)));
606233c2ff6SSatish Balay 
607233c2ff6SSatish Balay   PetscFunctionBegin;
60823ce1328SBarry Smith   if (fproc > size) fproc = size;
609e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
610e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
611233c2ff6SSatish Balay     else                         fproc++;
612233c2ff6SSatish Balay   }
613e005ede5SBarry Smith   *rank = fproc;
614233c2ff6SSatish Balay   PetscFunctionReturn(0);
615233c2ff6SSatish Balay }
616233c2ff6SSatish Balay #endif
617233c2ff6SSatish Balay 
618a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
619b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
6204a2ae208SSatish Balay #undef __FUNCT__
6214a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local"
622b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
623a2feddc5SSatish Balay {
624df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
625cf4f527aSSatish Balay   Mat            A = c->A;
626df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
6275d0c19d7SBarry Smith   const PetscInt **irow,**icol,*irow_i;
6285d0c19d7SBarry Smith   PetscInt       *nrow,*ncol,*w3,*w4,start;
6296849ba73SBarry Smith   PetscErrorCode ierr;
6309405f653SBarry Smith   PetscMPIInt    size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc;
6319405f653SBarry Smith   PetscInt       **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row;
632b24ad042SBarry Smith   PetscInt       nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol;
633b24ad042SBarry Smith   PetscInt       **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
6345d0c19d7SBarry Smith   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
635052f0c41SBarry Smith   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
636d0f46423SBarry Smith   PetscInt       bs=C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
637899cda47SBarry Smith   PetscInt       cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
638899cda47SBarry Smith   PetscInt       *bmap = c->garray,ctmp,rstart=c->rstartbs;
639d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
640d5b485abSSatish Balay   MPI_Request    *r_waits4,*s_waits3,*s_waits4;
641d5b485abSSatish Balay   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
642d5b485abSSatish Balay   MPI_Status     *r_status3,*r_status4,*s_status4;
643d5b485abSSatish Balay   MPI_Comm       comm;
6443eda8832SBarry Smith   MatScalar      **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB;
6453eda8832SBarry Smith   MatScalar      *a_a=a->a,*b_a=b->a;
646d4503c43SBarry Smith   PetscTruth     flag;
647b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1;
648c7dd2924SSatish Balay 
64980d608c0SSatish Balay #if defined (PETSC_USE_CTABLE)
650b24ad042SBarry Smith   PetscInt       tt;
651233c2ff6SSatish Balay   PetscTable     *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1;
652233c2ff6SSatish Balay #else
653b24ad042SBarry Smith   PetscInt       **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
654233c2ff6SSatish Balay #endif
655d5b485abSSatish Balay 
6563a40ed3dSBarry Smith   PetscFunctionBegin;
6577adad957SLisandro Dalcin   comm   = ((PetscObject)C)->comm;
6587adad957SLisandro Dalcin   tag0   = ((PetscObject)C)->tag;
659d5b485abSSatish Balay   size   = c->size;
660d5b485abSSatish Balay   rank   = c->rank;
661d5b485abSSatish Balay 
66234e6ae18SSatish Balay   /* Get some new tags to keep the communication clean */
66334e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
66434e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
66534e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
66634e6ae18SSatish Balay 
667052f0c41SBarry Smith #if defined(PETSC_USE_CTABLE)
66823bdbc58SBarry Smith   ierr = PetscMalloc4(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol);CHKERRQ(ierr);
669052f0c41SBarry Smith #else
67023bdbc58SBarry Smith   ierr = PetscMalloc5(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol,Mbs+1,PetscInt,&rtable);CHKERRQ(ierr);
671d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
672d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
673d5b485abSSatish Balay     jmax = c->rowners[i+1];
674d5b485abSSatish Balay     for (; j<jmax; j++) {
675d5b485abSSatish Balay       rtable[j] = i;
676d5b485abSSatish Balay     }
677d5b485abSSatish Balay   }
678233c2ff6SSatish Balay #endif
679233c2ff6SSatish Balay 
680233c2ff6SSatish Balay   for (i=0; i<ismax; i++) {
681233c2ff6SSatish Balay     ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
682233c2ff6SSatish Balay     ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
683233c2ff6SSatish Balay     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
684233c2ff6SSatish Balay     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
685233c2ff6SSatish Balay   }
686d5b485abSSatish Balay 
687d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg,and buffer space
688d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
68923bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr);
69023bdbc58SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
69123bdbc58SBarry Smith   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
69223bdbc58SBarry Smith   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
693d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
694b24ad042SBarry Smith     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
695d5b485abSSatish Balay     jmax   = nrow[i];
696d5b485abSSatish Balay     irow_i = irow[i];
697d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
698d5b485abSSatish Balay       row  = irow_i[j];
699233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
700899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
701233c2ff6SSatish Balay #else
702d5b485abSSatish Balay       proc = rtable[row];
703233c2ff6SSatish Balay #endif
704d5b485abSSatish Balay       w4[proc]++;
705d5b485abSSatish Balay     }
706d5b485abSSatish Balay     for (j=0; j<size; j++) {
707d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
708d5b485abSSatish Balay     }
709d5b485abSSatish Balay   }
710d5b485abSSatish Balay 
711d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
712e757655aSSatish Balay   msz      = 0;              /* total mesg length for all proc */
713d5b485abSSatish Balay   w1[rank] = 0;              /* no mesg sent to intself */
714d5b485abSSatish Balay   w3[rank] = 0;
715d5b485abSSatish Balay   for (i=0; i<size; i++) {
716d5b485abSSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
717d5b485abSSatish Balay   }
718b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/
719d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
720d5b485abSSatish Balay     if (w1[i]) { pa[j] = i; j++; }
721d5b485abSSatish Balay   }
722d5b485abSSatish Balay 
723d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
724d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
725d5b485abSSatish Balay     j     = pa[i];
726d5b485abSSatish Balay     w1[j] += w2[j] + 2* w3[j];
727d5b485abSSatish Balay     msz   += w1[j];
728d5b485abSSatish Balay   }
729d5b485abSSatish Balay 
730c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
731a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
732c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
733d5b485abSSatish Balay 
734c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
735c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
736c7dd2924SSatish Balay 
737c7dd2924SSatish Balay   ierr = PetscFree(onodes1);CHKERRQ(ierr);
738c7dd2924SSatish Balay   ierr = PetscFree(olengths1);CHKERRQ(ierr);
739d5b485abSSatish Balay 
740d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
74123bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscInt*,&sbuf1,size,PetscInt*,&ptr,2*msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr);
74223bdbc58SBarry Smith   ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
74323bdbc58SBarry Smith   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
744d5b485abSSatish Balay   {
745b24ad042SBarry Smith     PetscInt *iptr = tmp,ict = 0;
746d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
747d5b485abSSatish Balay       j         = pa[i];
748d5b485abSSatish Balay       iptr     += ict;
749d5b485abSSatish Balay       sbuf1[j]  = iptr;
750d5b485abSSatish Balay       ict       = w1[j];
751d5b485abSSatish Balay     }
752d5b485abSSatish Balay   }
753d5b485abSSatish Balay 
754d5b485abSSatish Balay   /* Form the outgoing messages */
755d5b485abSSatish Balay   /* Initialise the header space */
756d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
757d5b485abSSatish Balay     j           = pa[i];
758d5b485abSSatish Balay     sbuf1[j][0] = 0;
759b24ad042SBarry Smith     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
760d5b485abSSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
761d5b485abSSatish Balay   }
762d5b485abSSatish Balay 
763d5b485abSSatish Balay   /* Parse the isrow and copy data into outbuf */
764d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
765b24ad042SBarry Smith     ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
766d5b485abSSatish Balay     irow_i = irow[i];
767d5b485abSSatish Balay     jmax   = nrow[i];
768d5b485abSSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
769d5b485abSSatish Balay       row  = irow_i[j];
770233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
771899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
772233c2ff6SSatish Balay #else
773d5b485abSSatish Balay       proc = rtable[row];
774233c2ff6SSatish Balay #endif
775d5b485abSSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
776d5b485abSSatish Balay         ctr[proc]++;
777d5b485abSSatish Balay         *ptr[proc] = row;
778d5b485abSSatish Balay         ptr[proc]++;
779d5b485abSSatish Balay       }
780d5b485abSSatish Balay     }
781d5b485abSSatish Balay     /* Update the headers for the current IS */
782d5b485abSSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
78306ef35abSSatish Balay       if ((ctr_j = ctr[j])) {
784d5b485abSSatish Balay         sbuf1_j        = sbuf1[j];
785d5b485abSSatish Balay         k              = ++sbuf1_j[0];
786d5b485abSSatish Balay         sbuf1_j[2*k]   = ctr_j;
787d5b485abSSatish Balay         sbuf1_j[2*k-1] = i;
788d5b485abSSatish Balay       }
789d5b485abSSatish Balay     }
790d5b485abSSatish Balay   }
791d5b485abSSatish Balay 
792d5b485abSSatish Balay   /*  Now  post the sends */
793b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
794d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
795d5b485abSSatish Balay     j = pa[i];
796b24ad042SBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
797d5b485abSSatish Balay   }
798d5b485abSSatish Balay 
799d5b485abSSatish Balay   /* Post Recieves to capture the buffer size */
800b0a32e0cSBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr);
801b24ad042SBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr);
802d5b485abSSatish Balay   rbuf2[0] = tmp + msz;
803d5b485abSSatish Balay   for (i=1; i<nrqs; ++i) {
804d5b485abSSatish Balay     j        = pa[i];
805d5b485abSSatish Balay     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
806d5b485abSSatish Balay   }
807d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
808d5b485abSSatish Balay     j    = pa[i];
809b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
810d5b485abSSatish Balay   }
811d5b485abSSatish Balay 
812d5b485abSSatish Balay   /* Send to other procs the buf size they should allocate */
813d5b485abSSatish Balay 
814d5b485abSSatish Balay   /* Receive messages*/
815b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
816b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr);
81723bdbc58SBarry Smith   ierr = PetscMalloc3(nrqr+1,PetscInt*,&sbuf2,nrqr,PetscInt,&req_size,nrqr,PetscInt,&req_source);CHKERRQ(ierr);
818d5b485abSSatish Balay   {
819df36731bSBarry Smith     Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
820b24ad042SBarry Smith     PetscInt    *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
821d5b485abSSatish Balay 
822d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
823999d9058SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
824999d9058SBarry Smith       req_size[idex] = 0;
825999d9058SBarry Smith       rbuf1_i         = rbuf1[idex];
826d5b485abSSatish Balay       start           = 2*rbuf1_i[0] + 1;
827b24ad042SBarry Smith       ierr            = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
828b24ad042SBarry Smith       ierr            = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr);
829999d9058SBarry Smith       sbuf2_i         = sbuf2[idex];
830d5b485abSSatish Balay       for (j=start; j<end; j++) {
831d5b485abSSatish Balay         id               = rbuf1_i[j] - rstart;
832d5b485abSSatish Balay         ncols            = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
833d5b485abSSatish Balay         sbuf2_i[j]       = ncols;
834999d9058SBarry Smith         req_size[idex] += ncols;
835d5b485abSSatish Balay       }
836999d9058SBarry Smith       req_source[idex] = r_status1[i].MPI_SOURCE;
837d5b485abSSatish Balay       /* form the header */
838999d9058SBarry Smith       sbuf2_i[0]   = req_size[idex];
839d5b485abSSatish Balay       for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
840b24ad042SBarry Smith       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
841d5b485abSSatish Balay     }
842d5b485abSSatish Balay   }
843606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
844606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
845d5b485abSSatish Balay 
846d5b485abSSatish Balay   /*  recv buffer sizes */
847d5b485abSSatish Balay   /* Receive messages*/
848d5b485abSSatish Balay 
849b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr);
850b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr);
851b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr);
852b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr);
853b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr);
854d5b485abSSatish Balay 
855d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
856999d9058SBarry Smith     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
857b24ad042SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr);
858999d9058SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr);
859b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
860b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
861d5b485abSSatish Balay   }
862606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
863606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
864d5b485abSSatish Balay 
865d5b485abSSatish Balay   /* Wait on sends1 and sends2 */
866b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr);
867b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr);
868d5b485abSSatish Balay 
8690c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
8700c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
871606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
872606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
873606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
874606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
875d5b485abSSatish Balay 
876d5b485abSSatish Balay   /* Now allocate buffers for a->j, and send them off */
877b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr);
878d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
879b24ad042SBarry Smith   ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr);
880d5b485abSSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
881d5b485abSSatish Balay 
882b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr);
883d5b485abSSatish Balay   {
884d5b485abSSatish Balay      for (i=0; i<nrqr; i++) {
885d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
886d5b485abSSatish Balay       sbuf_aj_i = sbuf_aj[i];
887d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0] + 1;
888d5b485abSSatish Balay       ct2       = 0;
889d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
890d5b485abSSatish Balay         kmax = rbuf1[i][2*j];
891d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
892e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
893d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
894d5b485abSSatish Balay           ncols  = nzA + nzB;
895d5b485abSSatish Balay           cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
896d5b485abSSatish Balay 
897d5b485abSSatish Balay           /* load the column indices for this row into cols*/
898d5b485abSSatish Balay           cols  = sbuf_aj_i + ct2;
899d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
900d5b485abSSatish Balay             if ((ctmp = bmap[cworkB[l]]) < cstart)  cols[l] = ctmp;
901d5b485abSSatish Balay             else break;
902d5b485abSSatish Balay           }
903d5b485abSSatish Balay           imark = l;
904d5b485abSSatish Balay           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
905d5b485abSSatish Balay           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
906d5b485abSSatish Balay           ct2 += ncols;
907d5b485abSSatish Balay         }
908d5b485abSSatish Balay       }
909b24ad042SBarry Smith       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
910d5b485abSSatish Balay     }
911d5b485abSSatish Balay   }
912b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr);
913b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr);
914d5b485abSSatish Balay 
915d5b485abSSatish Balay   /* Allocate buffers for a->a, and send them off */
91682502324SSatish Balay   ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr);
917d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
91882502324SSatish Balay   ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr);
919a2feddc5SSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
920d5b485abSSatish Balay 
921b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr);
922d5b485abSSatish Balay   {
923d5b485abSSatish Balay     for (i=0; i<nrqr; i++) {
924d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
925d5b485abSSatish Balay       sbuf_aa_i = sbuf_aa[i];
926d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0]+1;
927d5b485abSSatish Balay       ct2       = 0;
928d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
929d5b485abSSatish Balay         kmax = rbuf1_i[2*j];
930d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
931e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
932d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
933d5b485abSSatish Balay           ncols  = nzA + nzB;
934d5b485abSSatish Balay           cworkA = a_j + a_i[row];     cworkB = b_j + b_i[row];
935a2feddc5SSatish Balay           vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2;
936d5b485abSSatish Balay 
937d5b485abSSatish Balay           /* load the column values for this row into vals*/
9385b83ace0SSatish Balay           vals  = sbuf_aa_i+ct2*bs2;
939d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
9403a40ed3dSBarry Smith             if ((bmap[cworkB[l]]) < cstart) {
9413eda8832SBarry Smith               ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9423a40ed3dSBarry Smith             }
943d5b485abSSatish Balay             else break;
944d5b485abSSatish Balay           }
945d5b485abSSatish Balay           imark = l;
9463a40ed3dSBarry Smith           for (l=0; l<nzA; l++) {
9473eda8832SBarry Smith             ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9483a40ed3dSBarry Smith           }
9493a40ed3dSBarry Smith           for (l=imark; l<nzB; l++) {
9503eda8832SBarry Smith             ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9513a40ed3dSBarry Smith           }
952d5b485abSSatish Balay           ct2 += ncols;
953d5b485abSSatish Balay         }
954d5b485abSSatish Balay       }
9553eda8832SBarry Smith       ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
956d5b485abSSatish Balay     }
957d5b485abSSatish Balay   }
958b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr);
959b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr);
960606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
961d5b485abSSatish Balay 
962d5b485abSSatish Balay   /* Form the matrix */
963d5b485abSSatish Balay   /* create col map */
964d5b485abSSatish Balay   {
9655d0c19d7SBarry Smith     const PetscInt *icol_i;
966233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
967233c2ff6SSatish Balay     /* Create row map*/
968b0a32e0cSBarry Smith     ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&colmaps);CHKERRQ(ierr);
969ff0794f7SSatish Balay     for (i=0; i<ismax; i++) {
970ff0794f7SSatish Balay       ierr = PetscTableCreate(ncol[i]+1,&colmaps[i]);CHKERRQ(ierr);
971233c2ff6SSatish Balay     }
972233c2ff6SSatish Balay #else
97323bdbc58SBarry Smith     ierr    = PetscMalloc(ismax*sizeof(PetscInt*),&cmap);CHKERRQ(ierr);
97423bdbc58SBarry Smith     ierr    = PetscMalloc(ismax*c->Nbs*sizeof(PetscInt),&cmap[0]);CHKERRQ(ierr);
975a2feddc5SSatish Balay     for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; }
976233c2ff6SSatish Balay #endif
977d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
978d5b485abSSatish Balay       jmax   = ncol[i];
979d5b485abSSatish Balay       icol_i = icol[i];
980233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
981233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
982233c2ff6SSatish Balay       for (j=0; j<jmax; j++) {
983001aedefSBarry Smith 	ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr);
984233c2ff6SSatish Balay       }
985233c2ff6SSatish Balay #else
986d5b485abSSatish Balay       cmap_i = cmap[i];
987d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
988d5b485abSSatish Balay         cmap_i[icol_i[j]] = j+1;
989d5b485abSSatish Balay       }
990233c2ff6SSatish Balay #endif
991d5b485abSSatish Balay     }
992d5b485abSSatish Balay   }
993d5b485abSSatish Balay 
994d5b485abSSatish Balay   /* Create lens which is required for MatCreate... */
995d5b485abSSatish Balay   for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
996052f0c41SBarry Smith   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr);
997b24ad042SBarry Smith   lens[0] = (PetscInt*)(lens + ismax);
998b24ad042SBarry Smith   ierr    = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr);
999d5b485abSSatish Balay   for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
1000d5b485abSSatish Balay 
1001d5b485abSSatish Balay   /* Update lens from local data */
1002d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1003d5b485abSSatish Balay     jmax   = nrow[i];
1004233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1005233c2ff6SSatish Balay     lcol1_gcol1 = colmaps[i];
1006233c2ff6SSatish Balay #else
1007d5b485abSSatish Balay     cmap_i = cmap[i];
1008233c2ff6SSatish Balay #endif
1009d5b485abSSatish Balay     irow_i = irow[i];
1010d5b485abSSatish Balay     lens_i = lens[i];
1011d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1012d5b485abSSatish Balay       row  = irow_i[j];
1013233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1014899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1015233c2ff6SSatish Balay #else
1016d5b485abSSatish Balay       proc = rtable[row];
1017233c2ff6SSatish Balay #endif
1018d5b485abSSatish Balay       if (proc == rank) {
101936f4e84dSSatish Balay         /* Get indices from matA and then from matB */
102036f4e84dSSatish Balay         row    = row - rstart;
102136f4e84dSSatish Balay         nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
102236f4e84dSSatish Balay         cworkA =  a_j + a_i[row]; cworkB = b_j + b_i[row];
1023233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1024233c2ff6SSatish Balay        for (k=0; k<nzA; k++) {
1025233c2ff6SSatish Balay 	 ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1026233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1027233c2ff6SSatish Balay         }
1028233c2ff6SSatish Balay         for (k=0; k<nzB; k++) {
1029233c2ff6SSatish Balay 	  ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1030233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1031233c2ff6SSatish Balay         }
1032233c2ff6SSatish Balay #else
1033ca161407SBarry Smith         for (k=0; k<nzA; k++) {
103436f4e84dSSatish Balay           if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; }
1035ca161407SBarry Smith         }
1036ca161407SBarry Smith         for (k=0; k<nzB; k++) {
103736f4e84dSSatish Balay           if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; }
1038d5b485abSSatish Balay         }
1039233c2ff6SSatish Balay #endif
1040d5b485abSSatish Balay       }
1041a2feddc5SSatish Balay     }
1042ca161407SBarry Smith   }
1043233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1044d5b485abSSatish Balay   /* Create row map*/
1045b0a32e0cSBarry Smith   ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rowmaps);CHKERRQ(ierr);
1046ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
1047ff0794f7SSatish Balay     ierr = PetscTableCreate(nrow[i]+1,&rowmaps[i]);CHKERRQ(ierr);
1048233c2ff6SSatish Balay   }
1049233c2ff6SSatish Balay #else
1050233c2ff6SSatish Balay   /* Create row map*/
1051052f0c41SBarry Smith   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr);
1052b24ad042SBarry Smith   rmap[0] = (PetscInt*)(rmap + ismax);
1053b24ad042SBarry Smith   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr);
1054233c2ff6SSatish Balay   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;}
1055233c2ff6SSatish Balay #endif
1056d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1057d5b485abSSatish Balay     irow_i = irow[i];
1058d5b485abSSatish Balay     jmax   = nrow[i];
1059233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1060233c2ff6SSatish Balay     lrow1_grow1 = rowmaps[i];
1061233c2ff6SSatish Balay     for (j=0; j<jmax; j++) {
1062233c2ff6SSatish Balay       ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr);
1063233c2ff6SSatish Balay     }
1064233c2ff6SSatish Balay #else
1065233c2ff6SSatish Balay     rmap_i = rmap[i];
1066d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1067d5b485abSSatish Balay       rmap_i[irow_i[j]] = j;
1068d5b485abSSatish Balay     }
1069233c2ff6SSatish Balay #endif
1070d5b485abSSatish Balay   }
1071d5b485abSSatish Balay 
1072d5b485abSSatish Balay   /* Update lens from offproc data */
1073d5b485abSSatish Balay   {
1074b24ad042SBarry Smith     PetscInt    *rbuf2_i,*rbuf3_i,*sbuf1_i;
1075b24ad042SBarry Smith     PetscMPIInt ii;
1076d5b485abSSatish Balay 
1077d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1078b24ad042SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr);
1079b24ad042SBarry Smith       idex   = pa[ii];
1080999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1081d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1082d5b485abSSatish Balay       ct1     = 2*jmax+1;
1083d5b485abSSatish Balay       ct2     = 0;
1084b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1085b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
1086d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1087d5b485abSSatish Balay         is_no   = sbuf1_i[2*j-1];
1088d5b485abSSatish Balay         max1    = sbuf1_i[2*j];
1089d5b485abSSatish Balay         lens_i  = lens[is_no];
1090233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1091233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1092233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1093233c2ff6SSatish Balay #else
1094d5b485abSSatish Balay         cmap_i  = cmap[is_no];
1095d5b485abSSatish Balay 	rmap_i  = rmap[is_no];
1096233c2ff6SSatish Balay #endif
1097d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1098233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1099233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1100233c2ff6SSatish Balay           row--;
1101e005ede5SBarry Smith           if (row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); }
1102233c2ff6SSatish Balay #else
1103d5b485abSSatish Balay           row  = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
1104233c2ff6SSatish Balay #endif
1105d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1106d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1107233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1108233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
1109233c2ff6SSatish Balay 	    if (tt) {
1110233c2ff6SSatish Balay               lens_i[row]++;
1111233c2ff6SSatish Balay             }
1112233c2ff6SSatish Balay #else
1113d5b485abSSatish Balay             if (cmap_i[rbuf3_i[ct2]]) {
1114d5b485abSSatish Balay               lens_i[row]++;
1115d5b485abSSatish Balay             }
1116233c2ff6SSatish Balay #endif
1117d5b485abSSatish Balay           }
1118d5b485abSSatish Balay         }
1119d5b485abSSatish Balay       }
1120d5b485abSSatish Balay     }
1121d5b485abSSatish Balay   }
1122606d414cSSatish Balay   ierr = PetscFree(r_status3);CHKERRQ(ierr);
1123606d414cSSatish Balay   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
11240c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
1125606d414cSSatish Balay   ierr = PetscFree(s_status3);CHKERRQ(ierr);
1126606d414cSSatish Balay   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1127d5b485abSSatish Balay 
1128d5b485abSSatish Balay   /* Create the submatrices */
1129d5b485abSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1130d5b485abSSatish Balay     /*
1131d5b485abSSatish Balay         Assumes new rows are same length as the old rows, hence bug!
1132d5b485abSSatish Balay     */
1133d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1134df36731bSBarry Smith       mat = (Mat_SeqBAIJ *)(submats[i]->data);
1135d0f46423SBarry Smith       if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) {
113629bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1137d5b485abSSatish Balay       }
1138b24ad042SBarry Smith       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
1139abc0a331SBarry Smith       if (!flag) {
114029bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1141d5b485abSSatish Balay       }
1142d5b485abSSatish Balay       /* Initial matrix as if empty */
1143b24ad042SBarry Smith       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr);
1144ce60de66SLois Curfman McInnes       submats[i]->factor = C->factor;
1145d5b485abSSatish Balay     }
1146ca161407SBarry Smith   } else {
1147d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1148f69a0ea3SMatthew Knepley       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
1149f69a0ea3SMatthew Knepley       ierr = MatSetSizes(submats[i],nrow[i]*bs,ncol[i]*bs,nrow[i]*bs,ncol[i]*bs);CHKERRQ(ierr);
11507adad957SLisandro Dalcin       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
1151d0f46423SBarry Smith       ierr = MatSeqBAIJSetPreallocation(submats[i],C->rmap->bs,0,lens[i]);CHKERRQ(ierr);
1152d0f46423SBarry Smith       ierr = MatSeqSBAIJSetPreallocation(submats[i],C->rmap->bs,0,lens[i]);CHKERRQ(ierr);
1153d5b485abSSatish Balay     }
1154d5b485abSSatish Balay   }
1155d5b485abSSatish Balay 
1156d5b485abSSatish Balay   /* Assemble the matrices */
1157d5b485abSSatish Balay   /* First assemble the local rows */
1158d5b485abSSatish Balay   {
1159b24ad042SBarry Smith     PetscInt       ilen_row,*imat_ilen,*imat_j,*imat_i;
11603eda8832SBarry Smith     MatScalar *imat_a;
1161d5b485abSSatish Balay 
1162d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1163df36731bSBarry Smith       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1164d5b485abSSatish Balay       imat_ilen = mat->ilen;
1165d5b485abSSatish Balay       imat_j    = mat->j;
1166d5b485abSSatish Balay       imat_i    = mat->i;
1167d5b485abSSatish Balay       imat_a    = mat->a;
1168233c2ff6SSatish Balay 
1169233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1170233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
1171233c2ff6SSatish Balay       lrow1_grow1 = rowmaps[i];
1172233c2ff6SSatish Balay #else
1173d5b485abSSatish Balay       cmap_i    = cmap[i];
1174d5b485abSSatish Balay       rmap_i    = rmap[i];
1175233c2ff6SSatish Balay #endif
1176d5b485abSSatish Balay       irow_i    = irow[i];
1177d5b485abSSatish Balay       jmax      = nrow[i];
1178d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1179d5b485abSSatish Balay         row      = irow_i[j];
1180233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1181899cda47SBarry Smith 	ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1182233c2ff6SSatish Balay #else
1183d5b485abSSatish Balay 	proc = rtable[row];
1184233c2ff6SSatish Balay #endif
1185d5b485abSSatish Balay         if (proc == rank) {
118636f4e84dSSatish Balay           row      = row - rstart;
118736f4e84dSSatish Balay           nzA      = a_i[row+1] - a_i[row];
118836f4e84dSSatish Balay           nzB      = b_i[row+1] - b_i[row];
118936f4e84dSSatish Balay           cworkA   = a_j + a_i[row];
119036f4e84dSSatish Balay           cworkB   = b_j + b_i[row];
119136f4e84dSSatish Balay           vworkA   = a_a + a_i[row]*bs2;
119236f4e84dSSatish Balay           vworkB   = b_a + b_i[row]*bs2;
1193233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1194233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr);
1195233c2ff6SSatish Balay           row--;
1196e005ede5SBarry Smith           if (row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); }
1197233c2ff6SSatish Balay #else
119836f4e84dSSatish Balay           row      = rmap_i[row + rstart];
1199233c2ff6SSatish Balay #endif
1200df36731bSBarry Smith           mat_i    = imat_i[row];
120136f4e84dSSatish Balay           mat_a    = imat_a + mat_i*bs2;
1202d5b485abSSatish Balay           mat_j    = imat_j + mat_i;
120336f4e84dSSatish Balay           ilen_row = imat_ilen[row];
120436f4e84dSSatish Balay 
120536f4e84dSSatish Balay           /* load the column indices for this row into cols*/
120636f4e84dSSatish Balay           for (l=0; l<nzB; l++) {
120736f4e84dSSatish Balay 	    if ((ctmp = bmap[cworkB[l]]) < cstart) {
1208233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1209233c2ff6SSatish Balay 	      ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr);
1210233c2ff6SSatish Balay 	      if (tcol) {
1211233c2ff6SSatish Balay #else
121236f4e84dSSatish Balay               if ((tcol = cmap_i[ctmp])) {
1213233c2ff6SSatish Balay #endif
1214df36731bSBarry Smith                 *mat_j++ = tcol - 1;
12153eda8832SBarry Smith                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1216549d3d68SSatish Balay                 mat_a   += bs2;
1217d5b485abSSatish Balay                 ilen_row++;
1218d5b485abSSatish Balay               }
1219ca161407SBarry Smith             } else break;
122036f4e84dSSatish Balay           }
122136f4e84dSSatish Balay           imark = l;
122236f4e84dSSatish Balay           for (l=0; l<nzA; l++) {
1223233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1224233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1225233c2ff6SSatish Balay 	    if (tcol) {
1226233c2ff6SSatish Balay #else
122736f4e84dSSatish Balay 	    if ((tcol = cmap_i[cstart + cworkA[l]])) {
1228233c2ff6SSatish Balay #endif
122936f4e84dSSatish Balay               *mat_j++ = tcol - 1;
12303eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1231549d3d68SSatish Balay               mat_a   += bs2;
123236f4e84dSSatish Balay               ilen_row++;
123336f4e84dSSatish Balay             }
123436f4e84dSSatish Balay           }
123536f4e84dSSatish Balay           for (l=imark; l<nzB; l++) {
1236233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1237233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1238233c2ff6SSatish Balay 	    if (tcol) {
1239233c2ff6SSatish Balay #else
124036f4e84dSSatish Balay             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1241233c2ff6SSatish Balay #endif
124236f4e84dSSatish Balay               *mat_j++ = tcol - 1;
12433eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1244549d3d68SSatish Balay               mat_a   += bs2;
124536f4e84dSSatish Balay               ilen_row++;
124636f4e84dSSatish Balay             }
124736f4e84dSSatish Balay           }
1248d5b485abSSatish Balay           imat_ilen[row] = ilen_row;
1249d5b485abSSatish Balay         }
1250d5b485abSSatish Balay       }
125136f4e84dSSatish Balay 
1252d5b485abSSatish Balay     }
1253d5b485abSSatish Balay   }
1254d5b485abSSatish Balay 
1255d5b485abSSatish Balay   /*   Now assemble the off proc rows*/
1256d5b485abSSatish Balay   {
1257b24ad042SBarry Smith     PetscInt    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1258b24ad042SBarry Smith     PetscInt    *imat_j,*imat_i;
12593eda8832SBarry Smith     MatScalar   *imat_a,*rbuf4_i;
1260b24ad042SBarry Smith     PetscMPIInt ii;
1261d5b485abSSatish Balay 
1262d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1263b24ad042SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr);
1264b24ad042SBarry Smith       idex   = pa[ii];
1265999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1266d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1267d5b485abSSatish Balay       ct1     = 2*jmax + 1;
1268d5b485abSSatish Balay       ct2     = 0;
1269b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1270b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
1271b24ad042SBarry Smith       rbuf4_i = rbuf4[ii];
1272d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1273d5b485abSSatish Balay         is_no     = sbuf1_i[2*j-1];
1274233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1275233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1276233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1277233c2ff6SSatish Balay #else
1278d5b485abSSatish Balay         rmap_i    = rmap[is_no];
1279d5b485abSSatish Balay         cmap_i    = cmap[is_no];
1280233c2ff6SSatish Balay #endif
1281df36731bSBarry Smith         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1282d5b485abSSatish Balay         imat_ilen = mat->ilen;
1283d5b485abSSatish Balay         imat_j    = mat->j;
1284d5b485abSSatish Balay         imat_i    = mat->i;
1285d5b485abSSatish Balay         imat_a    = mat->a;
1286d5b485abSSatish Balay         max1      = sbuf1_i[2*j];
1287d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1288d5b485abSSatish Balay           row   = sbuf1_i[ct1];
1289233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1290233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr);
1291233c2ff6SSatish Balay           row--;
1292e005ede5SBarry Smith           if(row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); }
1293233c2ff6SSatish Balay #else
1294d5b485abSSatish Balay           row   = rmap_i[row];
1295233c2ff6SSatish Balay #endif
1296d5b485abSSatish Balay           ilen  = imat_ilen[row];
1297df36731bSBarry Smith           mat_i = imat_i[row];
129836f4e84dSSatish Balay           mat_a = imat_a + mat_i*bs2;
1299d5b485abSSatish Balay           mat_j = imat_j + mat_i;
1300d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1301d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1302233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1303233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1304233c2ff6SSatish Balay 	    if (tcol) {
1305233c2ff6SSatish Balay #else
1306d5b485abSSatish Balay 	    if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1307233c2ff6SSatish Balay #endif
1308df36731bSBarry Smith               *mat_j++    = tcol - 1;
130936f4e84dSSatish Balay               /* *mat_a++= rbuf4_i[ct2]; */
13103eda8832SBarry Smith               ierr        = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1311549d3d68SSatish Balay               mat_a      += bs2;
1312d5b485abSSatish Balay               ilen++;
1313d5b485abSSatish Balay             }
1314d5b485abSSatish Balay           }
1315d5b485abSSatish Balay           imat_ilen[row] = ilen;
1316d5b485abSSatish Balay         }
1317d5b485abSSatish Balay       }
1318d5b485abSSatish Balay     }
1319d5b485abSSatish Balay   }
1320606d414cSSatish Balay   ierr = PetscFree(r_status4);CHKERRQ(ierr);
1321606d414cSSatish Balay   ierr = PetscFree(r_waits4);CHKERRQ(ierr);
13220c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
1323606d414cSSatish Balay   ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1324606d414cSSatish Balay   ierr = PetscFree(s_status4);CHKERRQ(ierr);
1325d5b485abSSatish Balay 
1326d5b485abSSatish Balay   /* Restore the indices */
1327d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1328d5b485abSSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1329d5b485abSSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1330d5b485abSSatish Balay   }
1331d5b485abSSatish Balay 
1332d5b485abSSatish Balay   /* Destroy allocated memory */
133323bdbc58SBarry Smith #if defined(PETSC_USE_CTABLE)
133423bdbc58SBarry Smith   ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr);
133523bdbc58SBarry Smith #else
133623bdbc58SBarry Smith   ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr);
133723bdbc58SBarry Smith #endif
133823bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
1339606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
1340d5b485abSSatish Balay 
134123bdbc58SBarry Smith   ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
1342606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1343606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1344d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
1345606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1346d5b485abSSatish Balay   }
1347d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1348606d414cSSatish Balay     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1349606d414cSSatish Balay     ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);
1350d5b485abSSatish Balay   }
135123bdbc58SBarry Smith   ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr);
1352606d414cSSatish Balay   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1353606d414cSSatish Balay   ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1354606d414cSSatish Balay   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1355606d414cSSatish Balay   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
1356606d414cSSatish Balay   ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1357606d414cSSatish Balay   ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
1358d5b485abSSatish Balay 
1359233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1360ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
13619c666560SBarry Smith     ierr = PetscTableDestroy(rowmaps[i]);CHKERRQ(ierr);
13629c666560SBarry Smith     ierr = PetscTableDestroy(colmaps[i]);CHKERRQ(ierr);
1363233c2ff6SSatish Balay   }
1364233c2ff6SSatish Balay   ierr = PetscFree(colmaps);CHKERRQ(ierr);
1365233c2ff6SSatish Balay   ierr = PetscFree(rowmaps);CHKERRQ(ierr);
1366233c2ff6SSatish Balay #else
1367606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
136823bdbc58SBarry Smith   ierr = PetscFree(cmap[0]);CHKERRQ(ierr);
1369233c2ff6SSatish Balay   ierr = PetscFree(cmap);CHKERRQ(ierr);
1370233c2ff6SSatish Balay #endif
1371606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1372d5b485abSSatish Balay 
1373d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
137436f4e84dSSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
137536f4e84dSSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1376d5b485abSSatish Balay   }
13773a40ed3dSBarry Smith   PetscFunctionReturn(0);
1378d5b485abSSatish Balay }
1379ca161407SBarry Smith 
1380