xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision db41cccff32efbb513e6f5a641add0fbc314552d)
123bdbc58SBarry Smith 
2d5b485abSSatish Balay /*
3d5b485abSSatish Balay    Routines to compute overlapping regions of a parallel MPI matrix
4d5b485abSSatish Balay   and to find submatrices that were shared across processors.
5d5b485abSSatish Balay */
6c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h>
7c6db04a5SJed Brown #include <petscbt.h>
8d5b485abSSatish Balay 
9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char **,PetscInt*,PetscInt**);
10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt **,PetscInt**,PetscInt*);
1109573ac7SBarry Smith extern PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
1209573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
13d5b485abSSatish Balay 
144a2ae208SSatish Balay #undef __FUNCT__
154a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ"
16b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov)
17d5b485abSSatish Balay {
186849ba73SBarry Smith   PetscErrorCode ierr;
19d0f46423SBarry Smith   PetscInt       i,N=C->cmap->N, bs=C->rmap->bs;
2036f4e84dSSatish Balay   IS             *is_new;
2136f4e84dSSatish Balay 
223a40ed3dSBarry Smith   PetscFunctionBegin;
2382502324SSatish Balay   ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr);
24df36731bSBarry Smith   /* Convert the indices into block format */
2527f478b1SHong Zhang   ierr = ISCompressIndicesGeneral(N,bs,imax,is,is_new);CHKERRQ(ierr);
26e32f2f54SBarry Smith   if (ov < 0){ SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
27d5b485abSSatish Balay   for (i=0; i<ov; ++i) {
2836f4e84dSSatish Balay     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr);
29d5b485abSSatish Balay   }
30ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);}
3127f478b1SHong Zhang   ierr = ISExpandIndicesGeneral(N,bs,imax,is_new,is);CHKERRQ(ierr);
32ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);}
33606d414cSSatish Balay   ierr = PetscFree(is_new);CHKERRQ(ierr);
343a40ed3dSBarry Smith   PetscFunctionReturn(0);
35d5b485abSSatish Balay }
36d5b485abSSatish Balay 
37d5b485abSSatish Balay /*
38d5b485abSSatish Balay   Sample message format:
39d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
400e9b0e7eSHong Zhang   to index sets is[1], is[5]
41d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
42d5b485abSSatish Balay   -----------
43d5b485abSSatish Balay   mesg [1] = 1 => is[1]
44d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
45d5b485abSSatish Balay   -----------
46d5b485abSSatish Balay   mesg [5] = 5  => is[5]
47d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
48d5b485abSSatish Balay   -----------
49d5b485abSSatish Balay   mesg [7]
500e9b0e7eSHong Zhang   mesg [n]  data(is[1])
51d5b485abSSatish Balay   -----------
52d5b485abSSatish Balay   mesg[n+1]
53d5b485abSSatish Balay   mesg[m]  data(is[5])
54d5b485abSSatish Balay   -----------
55d5b485abSSatish Balay 
56d5b485abSSatish Balay   Notes:
57d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
58d5b485abSSatish Balay   nrqr - no of requests recieved (which have to be or which have been processed
59d5b485abSSatish Balay */
604a2ae208SSatish Balay #undef __FUNCT__
614a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once"
62*db41cccfSHong Zhang PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[])
63d5b485abSSatish Balay {
64df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
655d0c19d7SBarry Smith   const PetscInt **idx,*idx_i;
6624c049a4SHong Zhang   PetscInt       *n,*w3,*w4,**data,len;
676849ba73SBarry Smith   PetscErrorCode ierr;
68b24ad042SBarry Smith   PetscMPIInt    size,rank,tag1,tag2,*w2,*w1,nrqr;
69245d216aSHong Zhang   PetscInt       Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr;
708f9f447aSHong Zhang   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p;
71b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1,*onodes2,*olengths2;
72f1af5d2fSBarry Smith   PetscBT        *table;
73d5b485abSSatish Balay   MPI_Comm       comm;
74d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
75d5b485abSSatish Balay   MPI_Status     *s_status,*recv_status;
768f9f447aSHong Zhang   char           *t_p;
77d5b485abSSatish Balay 
783a40ed3dSBarry Smith   PetscFunctionBegin;
797adad957SLisandro Dalcin   comm   = ((PetscObject)C)->comm;
80d5b485abSSatish Balay   size   = c->size;
81d5b485abSSatish Balay   rank   = c->rank;
82df36731bSBarry Smith   Mbs    = c->Mbs;
83d5b485abSSatish Balay 
84c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
85c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
86c7dd2924SSatish Balay 
8724c049a4SHong Zhang   ierr = PetscMalloc2(imax+1,const PetscInt*,&idx,imax,PetscInt,&n);CHKERRQ(ierr);
8824c049a4SHong Zhang 
89d5b485abSSatish Balay   for (i=0; i<imax; i++) {
90d5b485abSSatish Balay     ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr);
91b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
92d5b485abSSatish Balay   }
93d5b485abSSatish Balay 
94d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
95d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
9623bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr);
9723bdbc58SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
9823bdbc58SBarry Smith   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
9923bdbc58SBarry Smith   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
100d5b485abSSatish Balay   for (i=0; i<imax; i++) {
101b24ad042SBarry Smith     ierr  = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
102d5b485abSSatish Balay     idx_i = idx[i];
103d5b485abSSatish Balay     len   = n[i];
104d5b485abSSatish Balay     for (j=0; j<len; j++) {
105d5b485abSSatish Balay       row  = idx_i[j];
1066b41c64dSBarry Smith       if (row < 0) {
107e32f2f54SBarry Smith         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries");
1086b41c64dSBarry Smith       }
10924c049a4SHong Zhang       ierr = PetscLayoutFindOwner(C->rmap,row,&proc);CHKERRQ(ierr);
110d5b485abSSatish Balay       w4[proc]++;
111d5b485abSSatish Balay     }
112d5b485abSSatish Balay     for (j=0; j<size; j++){
113d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
114d5b485abSSatish Balay     }
115d5b485abSSatish Balay   }
116d5b485abSSatish Balay 
117d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
118d5b485abSSatish Balay   msz      = 0;              /* total mesg length (for all proc */
1190e9b0e7eSHong Zhang   w1[rank] = 0;              /* no mesg sent to itself */
120d5b485abSSatish Balay   w3[rank] = 0;
121d5b485abSSatish Balay   for (i=0; i<size; i++) {
122d5b485abSSatish Balay     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
123d5b485abSSatish Balay   }
124d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
125b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr);
126d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
127d5b485abSSatish Balay     if (w1[i]) {pa[j] = i; j++;}
128d5b485abSSatish Balay   }
129d5b485abSSatish Balay 
130d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
131d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
132d5b485abSSatish Balay     j      = pa[i];
133d5b485abSSatish Balay     w1[j] += w2[j] + 2*w3[j];
134d5b485abSSatish Balay     msz   += w1[j];
135d5b485abSSatish Balay   }
136d5b485abSSatish Balay 
137c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
138a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
139c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
140d5b485abSSatish Balay 
141c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
142c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr);
143d5b485abSSatish Balay 
144d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
14523bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscInt*,&outdat,size,PetscInt*,&ptr,msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr);
14623bdbc58SBarry Smith   ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr);
14723bdbc58SBarry Smith   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
148d5b485abSSatish Balay   {
149b24ad042SBarry Smith     PetscInt *iptr = tmp,ict  = 0;
150d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
151d5b485abSSatish Balay       j         = pa[i];
152d5b485abSSatish Balay       iptr     +=  ict;
153d5b485abSSatish Balay       outdat[j] = iptr;
154d5b485abSSatish Balay       ict       = w1[j];
155d5b485abSSatish Balay     }
156d5b485abSSatish Balay   }
157d5b485abSSatish Balay 
158d5b485abSSatish Balay   /* Form the outgoing messages */
159d5b485abSSatish Balay   /*plug in the headers*/
160d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
161d5b485abSSatish Balay     j            = pa[i];
162d5b485abSSatish Balay     outdat[j][0] = 0;
163b24ad042SBarry Smith     ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
164d5b485abSSatish Balay     ptr[j]       = outdat[j] + 2*w3[j] + 1;
165d5b485abSSatish Balay   }
166d5b485abSSatish Balay 
167d5b485abSSatish Balay   /* Memory for doing local proc's work*/
168d5b485abSSatish Balay   {
1698f9f447aSHong Zhang     ierr = PetscMalloc5(imax,PetscBT,&table, imax,PetscInt*,&data, imax,PetscInt,&isz,
1708f9f447aSHong Zhang                         Mbs*imax,PetscInt,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,char,&t_p);CHKERRQ(ierr);
1718f9f447aSHong Zhang     ierr = PetscMemzero(table,imax*sizeof(PetscBT));CHKERRQ(ierr);
1728f9f447aSHong Zhang     ierr = PetscMemzero(data,imax*sizeof(PetscInt*));CHKERRQ(ierr);
1738f9f447aSHong Zhang     ierr = PetscMemzero(isz,imax*sizeof(PetscInt));CHKERRQ(ierr);
1748f9f447aSHong Zhang     ierr = PetscMemzero(d_p,Mbs*imax*sizeof(PetscInt));CHKERRQ(ierr);
1758f9f447aSHong Zhang     ierr = PetscMemzero(t_p,(Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char));CHKERRQ(ierr);
176d5b485abSSatish Balay 
177bbd702dbSSatish Balay     for (i=0; i<imax; i++) {
178b6410449SSatish Balay       table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
179bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
180d5b485abSSatish Balay     }
181d5b485abSSatish Balay   }
182d5b485abSSatish Balay 
183d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
184d5b485abSSatish Balay   {
185b24ad042SBarry Smith     PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j;
186f1af5d2fSBarry Smith     PetscBT  table_i;
187d5b485abSSatish Balay 
188d5b485abSSatish Balay     for (i=0; i<imax; i++) {
189b24ad042SBarry Smith       ierr    = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
190d5b485abSSatish Balay       n_i     = n[i];
191d5b485abSSatish Balay       table_i = table[i];
192d5b485abSSatish Balay       idx_i   = idx[i];
193d5b485abSSatish Balay       data_i  = data[i];
194d5b485abSSatish Balay       isz_i   = isz[i];
195d5b485abSSatish Balay       for (j=0;  j<n_i; j++) {  /* parse the indices of each IS */
196d5b485abSSatish Balay         row  = idx_i[j];
19724c049a4SHong Zhang         ierr = PetscLayoutFindOwner(C->rmap,row,&proc);CHKERRQ(ierr);
198d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
199d5b485abSSatish Balay           ctr[proc]++;
200d5b485abSSatish Balay           *ptr[proc] = row;
201d5b485abSSatish Balay           ptr[proc]++;
202d6b45a43SBarry Smith         } else { /* Update the local table */
203f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
204d5b485abSSatish Balay         }
205d5b485abSSatish Balay       }
206d5b485abSSatish Balay       /* Update the headers for the current IS */
207d5b485abSSatish Balay       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
208d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
209d5b485abSSatish Balay           outdat_j        = outdat[j];
210d5b485abSSatish Balay           k               = ++outdat_j[0];
211d5b485abSSatish Balay           outdat_j[2*k]   = ctr_j;
212d5b485abSSatish Balay           outdat_j[2*k-1] = i;
213d5b485abSSatish Balay         }
214d5b485abSSatish Balay       }
215d5b485abSSatish Balay       isz[i] = isz_i;
216d5b485abSSatish Balay     }
217d5b485abSSatish Balay   }
218d5b485abSSatish Balay 
219d5b485abSSatish Balay   /*  Now  post the sends */
220b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
221d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
222d5b485abSSatish Balay     j    = pa[i];
223b24ad042SBarry Smith     ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr);
224d5b485abSSatish Balay   }
225d5b485abSSatish Balay 
226d5b485abSSatish Balay   /* No longer need the original indices*/
227d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
228d5b485abSSatish Balay     ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr);
229d5b485abSSatish Balay   }
23024c049a4SHong Zhang   ierr = PetscFree2(idx,n);CHKERRQ(ierr);
231d5b485abSSatish Balay 
232d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
233d5b485abSSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
234d5b485abSSatish Balay   }
235d5b485abSSatish Balay 
236d5b485abSSatish Balay   /* Do Local work*/
237df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
238d5b485abSSatish Balay 
239d5b485abSSatish Balay   /* Receive messages*/
240b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr);
2410c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);}
242d5b485abSSatish Balay 
243b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr);
2440c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);}
245d5b485abSSatish Balay 
246d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
24723bdbc58SBarry Smith   ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr);
24823bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
249d5b485abSSatish Balay 
250b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&xdata);CHKERRQ(ierr);
251b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt),&isz1);CHKERRQ(ierr);
252df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
253c05d87d6SBarry Smith   ierr = PetscFree(rbuf[0]);CHKERRQ(ierr);
254606d414cSSatish Balay   ierr = PetscFree(rbuf);CHKERRQ(ierr);
255d5b485abSSatish Balay 
256d5b485abSSatish Balay   /* Send the data back*/
257d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
258d5b485abSSatish Balay   {
259b24ad042SBarry Smith     PetscMPIInt *rw1;
260d5b485abSSatish Balay 
261b24ad042SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&rw1);CHKERRQ(ierr);
262b24ad042SBarry Smith     ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr);
263c7dd2924SSatish Balay 
264d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
265d5b485abSSatish Balay       proc      = recv_status[i].MPI_SOURCE;
266e32f2f54SBarry Smith       if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch");
267d5b485abSSatish Balay       rw1[proc] = isz1[i];
268d5b485abSSatish Balay     }
269d5b485abSSatish Balay 
270c7dd2924SSatish Balay     ierr = PetscFree(onodes1);CHKERRQ(ierr);
271c7dd2924SSatish Balay     ierr = PetscFree(olengths1);CHKERRQ(ierr);
272d5b485abSSatish Balay 
273c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
274c7dd2924SSatish Balay     ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr);
27503512dc5SSatish Balay     ierr = PetscFree(rw1);CHKERRQ(ierr);
276c7dd2924SSatish Balay   }
277c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
278c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr);
279d5b485abSSatish Balay 
280d5b485abSSatish Balay   /*  Now  post the sends */
281b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
282d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
283d5b485abSSatish Balay     j    = recv_status[i].MPI_SOURCE;
284b24ad042SBarry Smith     ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr);
285d5b485abSSatish Balay   }
286d5b485abSSatish Balay 
287d5b485abSSatish Balay   /* receive work done on other processors*/
288d5b485abSSatish Balay   {
289b24ad042SBarry Smith     PetscMPIInt idex;
290b24ad042SBarry Smith     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
291f1af5d2fSBarry Smith     PetscBT     table_i;
292d5b485abSSatish Balay     MPI_Status  *status2;
293d5b485abSSatish Balay 
294169449a0SSatish Balay     ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr);
295d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
296999d9058SBarry Smith       ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr);
297d5b485abSSatish Balay       /* Process the message*/
298999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
299d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
300999d9058SBarry Smith       jmax    = rbuf2[idex][0];
301d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
302d5b485abSSatish Balay         max     = rbuf2_i[2*j];
303d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
304d5b485abSSatish Balay         isz_i   = isz[is_no];
305d5b485abSSatish Balay         data_i  = data[is_no];
306d5b485abSSatish Balay         table_i = table[is_no];
307d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
308d5b485abSSatish Balay           row = rbuf2_i[ct1];
309f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
310d5b485abSSatish Balay         }
311d5b485abSSatish Balay         isz[is_no] = isz_i;
312d5b485abSSatish Balay       }
313d5b485abSSatish Balay     }
3140c468ba9SBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);}
315606d414cSSatish Balay     ierr = PetscFree(status2);CHKERRQ(ierr);
316d5b485abSSatish Balay   }
317d5b485abSSatish Balay 
318d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
31970b3c8c7SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
320d5b485abSSatish Balay   }
321d5b485abSSatish Balay 
322c7dd2924SSatish Balay 
323c7dd2924SSatish Balay   ierr = PetscFree(onodes2);CHKERRQ(ierr);
324c7dd2924SSatish Balay   ierr = PetscFree(olengths2);CHKERRQ(ierr);
325c7dd2924SSatish Balay 
326606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
327c05d87d6SBarry Smith   ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr);
328606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
329606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
330606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
331606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
332606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
3338f9f447aSHong Zhang   ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr);
334606d414cSSatish Balay   ierr = PetscFree(s_status);CHKERRQ(ierr);
335606d414cSSatish Balay   ierr = PetscFree(recv_status);CHKERRQ(ierr);
336606d414cSSatish Balay   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
337606d414cSSatish Balay   ierr = PetscFree(xdata);CHKERRQ(ierr);
338606d414cSSatish Balay   ierr = PetscFree(isz1);CHKERRQ(ierr);
3393a40ed3dSBarry Smith   PetscFunctionReturn(0);
340d5b485abSSatish Balay }
341d5b485abSSatish Balay 
3424a2ae208SSatish Balay #undef __FUNCT__
3434a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local"
344d5b485abSSatish Balay /*
345df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
346d5b485abSSatish Balay        the work on the local processor.
347d5b485abSSatish Balay 
348d5b485abSSatish Balay      Inputs:
349df36731bSBarry Smith       C      - MAT_MPIBAIJ;
350d5b485abSSatish Balay       imax - total no of index sets processed at a time;
351df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
352d5b485abSSatish Balay 
353d5b485abSSatish Balay      Output:
3540e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
355d5b485abSSatish Balay                to each index set;
356d5b485abSSatish Balay       data   - pointer to the solutions
357d5b485abSSatish Balay */
358b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data)
359d5b485abSSatish Balay {
360df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
361d5b485abSSatish Balay   Mat         A = c->A,B = c->B;
362df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
363b24ad042SBarry Smith   PetscInt    start,end,val,max,rstart,cstart,*ai,*aj;
364b24ad042SBarry Smith   PetscInt    *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
365f1af5d2fSBarry Smith   PetscBT     table_i;
366d5b485abSSatish Balay 
3673a40ed3dSBarry Smith   PetscFunctionBegin;
368899cda47SBarry Smith   rstart = c->rstartbs;
369899cda47SBarry Smith   cstart = c->cstartbs;
370d5b485abSSatish Balay   ai     = a->i;
371df36731bSBarry Smith   aj     = a->j;
372d5b485abSSatish Balay   bi     = b->i;
373df36731bSBarry Smith   bj     = b->j;
374d5b485abSSatish Balay   garray = c->garray;
375d5b485abSSatish Balay 
376d5b485abSSatish Balay 
377d5b485abSSatish Balay   for (i=0; i<imax; i++) {
378d5b485abSSatish Balay     data_i  = data[i];
379d5b485abSSatish Balay     table_i = table[i];
380d5b485abSSatish Balay     isz_i   = isz[i];
381d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
382d5b485abSSatish Balay       row   = data_i[j] - rstart;
383d5b485abSSatish Balay       start = ai[row];
384d5b485abSSatish Balay       end   = ai[row+1];
385d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
386df36731bSBarry Smith         val = aj[k] + cstart;
387f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
388d5b485abSSatish Balay       }
389d5b485abSSatish Balay       start = bi[row];
390d5b485abSSatish Balay       end   = bi[row+1];
391d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
392df36731bSBarry Smith         val = garray[bj[k]];
393f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
394d5b485abSSatish Balay       }
395d5b485abSSatish Balay     }
396d5b485abSSatish Balay     isz[i] = isz_i;
397d5b485abSSatish Balay   }
3983a40ed3dSBarry Smith   PetscFunctionReturn(0);
399d5b485abSSatish Balay }
4004a2ae208SSatish Balay #undef __FUNCT__
4014a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive"
402d5b485abSSatish Balay /*
403df36731bSBarry Smith       MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages,
404d5b485abSSatish Balay          and return the output
405d5b485abSSatish Balay 
406d5b485abSSatish Balay          Input:
407d5b485abSSatish Balay            C    - the matrix
408d5b485abSSatish Balay            nrqr - no of messages being processed.
409d5b485abSSatish Balay            rbuf - an array of pointers to the recieved requests
410d5b485abSSatish Balay 
411d5b485abSSatish Balay          Output:
412d5b485abSSatish Balay            xdata - array of messages to be sent back
413d5b485abSSatish Balay            isz1  - size of each message
414d5b485abSSatish Balay 
415a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
416d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
4170e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of
418d5b485abSSatish Balay memory is used.
419d5b485abSSatish Balay 
420d5b485abSSatish Balay */
421b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
422d5b485abSSatish Balay {
423df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
424d5b485abSSatish Balay   Mat            A = c->A,B = c->B;
425df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
4266849ba73SBarry Smith   PetscErrorCode ierr;
427b24ad042SBarry Smith   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
428b24ad042SBarry Smith   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
429d2a212eaSBarry Smith   PetscInt       val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
430b24ad042SBarry Smith   PetscInt       *rbuf_i,kmax,rbuf_0;
431f1af5d2fSBarry Smith   PetscBT        xtable;
432d5b485abSSatish Balay 
4333a40ed3dSBarry Smith   PetscFunctionBegin;
434df36731bSBarry Smith   Mbs    = c->Mbs;
435899cda47SBarry Smith   rstart = c->rstartbs;
436899cda47SBarry Smith   cstart = c->cstartbs;
437d5b485abSSatish Balay   ai     = a->i;
438df36731bSBarry Smith   aj     = a->j;
439d5b485abSSatish Balay   bi     = b->i;
440df36731bSBarry Smith   bj     = b->j;
441d5b485abSSatish Balay   garray = c->garray;
442d5b485abSSatish Balay 
443d5b485abSSatish Balay 
444d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
445d5b485abSSatish Balay     rbuf_i  =  rbuf[i];
446d5b485abSSatish Balay     rbuf_0  =  rbuf_i[0];
447d5b485abSSatish Balay     ct     += rbuf_0;
448d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; }
449d5b485abSSatish Balay   }
450d5b485abSSatish Balay 
451701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
452701b8814SBarry Smith   else        max1 = 1;
453d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
454b24ad042SBarry Smith   ierr         = PetscMalloc(mem_estimate*sizeof(PetscInt),&xdata[0]);CHKERRQ(ierr);
455d5b485abSSatish Balay   ++no_malloc;
4566831982aSBarry Smith   ierr         = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr);
457b24ad042SBarry Smith   ierr         = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr);
458d5b485abSSatish Balay 
459d5b485abSSatish Balay   ct3 = 0;
460d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
461d5b485abSSatish Balay     rbuf_i =  rbuf[i];
462d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
463d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
464d5b485abSSatish Balay     ct2    =  ct1;
465d5b485abSSatish Balay     ct3    += ct1;
466d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
4676831982aSBarry Smith       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
468d5b485abSSatish Balay       oct2 = ct2;
469d5b485abSSatish Balay       kmax = rbuf_i[2*j];
470d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
471d5b485abSSatish Balay         row = rbuf_i[ct1];
472f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
473d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
474b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
475b24ad042SBarry Smith             ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
476b24ad042SBarry Smith             ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
477606d414cSSatish Balay             ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
478d5b485abSSatish Balay             xdata[0]     = tmp;
479d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
480d5b485abSSatish Balay             for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
481d5b485abSSatish Balay           }
482d5b485abSSatish Balay           xdata[i][ct2++] = row;
483d5b485abSSatish Balay           ct3++;
484d5b485abSSatish Balay         }
485d5b485abSSatish Balay       }
486d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
487d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
488d5b485abSSatish Balay         start = ai[row];
489d5b485abSSatish Balay         end   = ai[row+1];
490d5b485abSSatish Balay         for (l=start; l<end; l++) {
491df36731bSBarry Smith           val = aj[l] + cstart;
492f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
493d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
494b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
495b24ad042SBarry Smith               ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
496b24ad042SBarry Smith               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
497606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
498d5b485abSSatish Balay               xdata[0]     = tmp;
499d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
500d5b485abSSatish Balay               for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
501d5b485abSSatish Balay             }
502d5b485abSSatish Balay             xdata[i][ct2++] = val;
503d5b485abSSatish Balay             ct3++;
504d5b485abSSatish Balay           }
505d5b485abSSatish Balay         }
506d5b485abSSatish Balay         start = bi[row];
507d5b485abSSatish Balay         end   = bi[row+1];
508d5b485abSSatish Balay         for (l=start; l<end; l++) {
509df36731bSBarry Smith           val = garray[bj[l]];
510f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
511d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
512b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
513b24ad042SBarry Smith               ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr);
514b24ad042SBarry Smith               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
515606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
516d5b485abSSatish Balay               xdata[0]     = tmp;
517d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
518d5b485abSSatish Balay               for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
519d5b485abSSatish Balay             }
520d5b485abSSatish Balay             xdata[i][ct2++] = val;
521d5b485abSSatish Balay             ct3++;
522d5b485abSSatish Balay           }
523d5b485abSSatish Balay         }
524d5b485abSSatish Balay       }
525d5b485abSSatish Balay       /* Update the header*/
526d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
527d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
528d5b485abSSatish Balay     }
529d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
530d5b485abSSatish Balay     xdata[i+1]  = xdata[i] + ct2;
531d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
532d5b485abSSatish Balay   }
5336831982aSBarry Smith   ierr = PetscBTDestroy(xtable);CHKERRQ(ierr);
5341e2582c4SBarry Smith   ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr);
5353a40ed3dSBarry Smith   PetscFunctionReturn(0);
536d5b485abSSatish Balay }
537d5b485abSSatish Balay 
538b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *);
539a2feddc5SSatish Balay 
5404a2ae208SSatish Balay #undef __FUNCT__
5414a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ"
542b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
543d5b485abSSatish Balay {
54436f4e84dSSatish Balay   IS             *isrow_new,*iscol_new;
545cf4f527aSSatish Balay   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
5466849ba73SBarry Smith   PetscErrorCode ierr;
547d0f46423SBarry Smith   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,N=C->cmap->N,bs=C->rmap->bs;
548a2feddc5SSatish Balay 
5493a40ed3dSBarry Smith   PetscFunctionBegin;
5503a40ed3dSBarry Smith   /* The compression and expansion should be avoided. Does'nt point
551a2feddc5SSatish Balay      out errors might change the indices hence buggey */
552a2feddc5SSatish Balay 
55323bdbc58SBarry Smith   ierr = PetscMalloc2(ismax+1,IS,&isrow_new,ismax+1,IS,&iscol_new);CHKERRQ(ierr);
554d4503c43SBarry Smith   ierr = ISCompressIndicesGeneral(N,bs,ismax,isrow,isrow_new);CHKERRQ(ierr);
555d4503c43SBarry Smith   ierr = ISCompressIndicesGeneral(N,bs,ismax,iscol,iscol_new);CHKERRQ(ierr);
556cf4f527aSSatish Balay 
557cf4f527aSSatish Balay   /* Allocate memory to hold all the submatrices */
558cf4f527aSSatish Balay   if (scall != MAT_REUSE_MATRIX) {
55982502324SSatish Balay     ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr);
560cf4f527aSSatish Balay   }
561cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
562b24ad042SBarry Smith   nmax          = 20*1000000 / (c->Nbs * sizeof(PetscInt));
563329f5518SBarry Smith   if (!nmax) nmax = 1;
564cf4f527aSSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
565cf4f527aSSatish Balay 
566653e4784SBarry Smith   /* Make sure every processor loops through the nstages */
5677adad957SLisandro Dalcin   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);CHKERRQ(ierr);
568cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
569cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
570cf4f527aSSatish Balay     else if (pos == ismax) max_no = 0;
571cf4f527aSSatish Balay     else                   max_no = ismax-pos;
572cf4f527aSSatish Balay     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr);
573cf4f527aSSatish Balay     pos += max_no;
574cf4f527aSSatish Balay   }
57536f4e84dSSatish Balay 
57636f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
577ca161407SBarry Smith     ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr);
578ca161407SBarry Smith     ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr);
57936f4e84dSSatish Balay   }
58023bdbc58SBarry Smith   ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr);
5813a40ed3dSBarry Smith   PetscFunctionReturn(0);
582a2feddc5SSatish Balay }
583a2feddc5SSatish Balay 
584233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
5854a2ae208SSatish Balay #undef __FUNCT__
5864a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc"
587e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
588233c2ff6SSatish Balay {
589e005ede5SBarry Smith   PetscInt    nGlobalNd = proc_gnode[size];
590b9f7ace7SBarry Smith   PetscMPIInt fproc = PetscMPIIntCast( (PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)));
591233c2ff6SSatish Balay 
592233c2ff6SSatish Balay   PetscFunctionBegin;
59323ce1328SBarry Smith   if (fproc > size) fproc = size;
594e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
595e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
596233c2ff6SSatish Balay     else                         fproc++;
597233c2ff6SSatish Balay   }
598e005ede5SBarry Smith   *rank = fproc;
599233c2ff6SSatish Balay   PetscFunctionReturn(0);
600233c2ff6SSatish Balay }
601233c2ff6SSatish Balay #endif
602233c2ff6SSatish Balay 
603a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
604b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
6054a2ae208SSatish Balay #undef __FUNCT__
6064a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local"
607b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
608a2feddc5SSatish Balay {
609df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
610cf4f527aSSatish Balay   Mat            A = c->A;
611df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
6125d0c19d7SBarry Smith   const PetscInt **irow,**icol,*irow_i;
6135d0c19d7SBarry Smith   PetscInt       *nrow,*ncol,*w3,*w4,start;
6146849ba73SBarry Smith   PetscErrorCode ierr;
6159405f653SBarry Smith   PetscMPIInt    size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc;
6169405f653SBarry Smith   PetscInt       **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row;
617b24ad042SBarry Smith   PetscInt       nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol;
618b24ad042SBarry Smith   PetscInt       **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
6195d0c19d7SBarry Smith   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
620052f0c41SBarry Smith   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
621d0f46423SBarry Smith   PetscInt       bs=C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
622899cda47SBarry Smith   PetscInt       cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
623899cda47SBarry Smith   PetscInt       *bmap = c->garray,ctmp,rstart=c->rstartbs;
624d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
625d5b485abSSatish Balay   MPI_Request    *r_waits4,*s_waits3,*s_waits4;
626d5b485abSSatish Balay   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
627d5b485abSSatish Balay   MPI_Status     *r_status3,*r_status4,*s_status4;
628d5b485abSSatish Balay   MPI_Comm       comm;
6293eda8832SBarry Smith   MatScalar      **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB;
6303eda8832SBarry Smith   MatScalar      *a_a=a->a,*b_a=b->a;
631ace3abfcSBarry Smith   PetscBool      flag;
632b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1;
633c7dd2924SSatish Balay 
63480d608c0SSatish Balay #if defined (PETSC_USE_CTABLE)
635b24ad042SBarry Smith   PetscInt       tt;
636233c2ff6SSatish Balay   PetscTable     *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1;
637233c2ff6SSatish Balay #else
638b24ad042SBarry Smith   PetscInt       **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
639233c2ff6SSatish Balay #endif
640d5b485abSSatish Balay 
6413a40ed3dSBarry Smith   PetscFunctionBegin;
6427adad957SLisandro Dalcin   comm   = ((PetscObject)C)->comm;
6437adad957SLisandro Dalcin   tag0   = ((PetscObject)C)->tag;
644d5b485abSSatish Balay   size   = c->size;
645d5b485abSSatish Balay   rank   = c->rank;
646d5b485abSSatish Balay 
64734e6ae18SSatish Balay   /* Get some new tags to keep the communication clean */
64834e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
64934e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
65034e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
65134e6ae18SSatish Balay 
652052f0c41SBarry Smith #if defined(PETSC_USE_CTABLE)
65323bdbc58SBarry Smith   ierr = PetscMalloc4(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol);CHKERRQ(ierr);
654052f0c41SBarry Smith #else
65523bdbc58SBarry Smith   ierr = PetscMalloc5(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol,Mbs+1,PetscInt,&rtable);CHKERRQ(ierr);
656d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
657d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
658d5b485abSSatish Balay     jmax = c->rowners[i+1];
659d5b485abSSatish Balay     for (; j<jmax; j++) {
660d5b485abSSatish Balay       rtable[j] = i;
661d5b485abSSatish Balay     }
662d5b485abSSatish Balay   }
663233c2ff6SSatish Balay #endif
664233c2ff6SSatish Balay 
665233c2ff6SSatish Balay   for (i=0; i<ismax; i++) {
666233c2ff6SSatish Balay     ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
667233c2ff6SSatish Balay     ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
668233c2ff6SSatish Balay     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
669233c2ff6SSatish Balay     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
670233c2ff6SSatish Balay   }
671d5b485abSSatish Balay 
672d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg,and buffer space
673d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
67423bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr);
67523bdbc58SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
67623bdbc58SBarry Smith   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
67723bdbc58SBarry Smith   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
678d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
679b24ad042SBarry Smith     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
680d5b485abSSatish Balay     jmax   = nrow[i];
681d5b485abSSatish Balay     irow_i = irow[i];
682d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
683d5b485abSSatish Balay       row  = irow_i[j];
684233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
685899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
686233c2ff6SSatish Balay #else
687d5b485abSSatish Balay       proc = rtable[row];
688233c2ff6SSatish Balay #endif
689d5b485abSSatish Balay       w4[proc]++;
690d5b485abSSatish Balay     }
691d5b485abSSatish Balay     for (j=0; j<size; j++) {
692d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
693d5b485abSSatish Balay     }
694d5b485abSSatish Balay   }
695d5b485abSSatish Balay 
696d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
697e757655aSSatish Balay   msz      = 0;              /* total mesg length for all proc */
698d5b485abSSatish Balay   w1[rank] = 0;              /* no mesg sent to intself */
699d5b485abSSatish Balay   w3[rank] = 0;
700d5b485abSSatish Balay   for (i=0; i<size; i++) {
701d5b485abSSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
702d5b485abSSatish Balay   }
703b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/
704d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
705d5b485abSSatish Balay     if (w1[i]) { pa[j] = i; j++; }
706d5b485abSSatish Balay   }
707d5b485abSSatish Balay 
708d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
709d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
710d5b485abSSatish Balay     j     = pa[i];
711d5b485abSSatish Balay     w1[j] += w2[j] + 2* w3[j];
712d5b485abSSatish Balay     msz   += w1[j];
713d5b485abSSatish Balay   }
714d5b485abSSatish Balay 
715c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
716a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
717c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
718d5b485abSSatish Balay 
719c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
720c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
721c7dd2924SSatish Balay 
722c7dd2924SSatish Balay   ierr = PetscFree(onodes1);CHKERRQ(ierr);
723c7dd2924SSatish Balay   ierr = PetscFree(olengths1);CHKERRQ(ierr);
724d5b485abSSatish Balay 
725d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
72623bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscInt*,&sbuf1,size,PetscInt*,&ptr,2*msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr);
72723bdbc58SBarry Smith   ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
72823bdbc58SBarry Smith   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
729d5b485abSSatish Balay   {
730b24ad042SBarry Smith     PetscInt *iptr = tmp,ict = 0;
731d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
732d5b485abSSatish Balay       j         = pa[i];
733d5b485abSSatish Balay       iptr     += ict;
734d5b485abSSatish Balay       sbuf1[j]  = iptr;
735d5b485abSSatish Balay       ict       = w1[j];
736d5b485abSSatish Balay     }
737d5b485abSSatish Balay   }
738d5b485abSSatish Balay 
739d5b485abSSatish Balay   /* Form the outgoing messages */
740d5b485abSSatish Balay   /* Initialise the header space */
741d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
742d5b485abSSatish Balay     j           = pa[i];
743d5b485abSSatish Balay     sbuf1[j][0] = 0;
744b24ad042SBarry Smith     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
745d5b485abSSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
746d5b485abSSatish Balay   }
747d5b485abSSatish Balay 
748d5b485abSSatish Balay   /* Parse the isrow and copy data into outbuf */
749d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
750b24ad042SBarry Smith     ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
751d5b485abSSatish Balay     irow_i = irow[i];
752d5b485abSSatish Balay     jmax   = nrow[i];
753d5b485abSSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
754d5b485abSSatish Balay       row  = irow_i[j];
755233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
756899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
757233c2ff6SSatish Balay #else
758d5b485abSSatish Balay       proc = rtable[row];
759233c2ff6SSatish Balay #endif
760d5b485abSSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
761d5b485abSSatish Balay         ctr[proc]++;
762d5b485abSSatish Balay         *ptr[proc] = row;
763d5b485abSSatish Balay         ptr[proc]++;
764d5b485abSSatish Balay       }
765d5b485abSSatish Balay     }
766d5b485abSSatish Balay     /* Update the headers for the current IS */
767d5b485abSSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
76806ef35abSSatish Balay       if ((ctr_j = ctr[j])) {
769d5b485abSSatish Balay         sbuf1_j        = sbuf1[j];
770d5b485abSSatish Balay         k              = ++sbuf1_j[0];
771d5b485abSSatish Balay         sbuf1_j[2*k]   = ctr_j;
772d5b485abSSatish Balay         sbuf1_j[2*k-1] = i;
773d5b485abSSatish Balay       }
774d5b485abSSatish Balay     }
775d5b485abSSatish Balay   }
776d5b485abSSatish Balay 
777d5b485abSSatish Balay   /*  Now  post the sends */
778b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
779d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
780d5b485abSSatish Balay     j = pa[i];
781b24ad042SBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
782d5b485abSSatish Balay   }
783d5b485abSSatish Balay 
784d5b485abSSatish Balay   /* Post Recieves to capture the buffer size */
785b0a32e0cSBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr);
786b24ad042SBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr);
787d5b485abSSatish Balay   rbuf2[0] = tmp + msz;
788d5b485abSSatish Balay   for (i=1; i<nrqs; ++i) {
789d5b485abSSatish Balay     j        = pa[i];
790d5b485abSSatish Balay     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
791d5b485abSSatish Balay   }
792d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
793d5b485abSSatish Balay     j    = pa[i];
794b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
795d5b485abSSatish Balay   }
796d5b485abSSatish Balay 
797d5b485abSSatish Balay   /* Send to other procs the buf size they should allocate */
798d5b485abSSatish Balay 
799d5b485abSSatish Balay   /* Receive messages*/
800b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
801b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr);
80223bdbc58SBarry Smith   ierr = PetscMalloc3(nrqr+1,PetscInt*,&sbuf2,nrqr,PetscInt,&req_size,nrqr,PetscInt,&req_source);CHKERRQ(ierr);
803d5b485abSSatish Balay   {
804df36731bSBarry Smith     Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
805b24ad042SBarry Smith     PetscInt    *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
806d5b485abSSatish Balay 
807d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
808999d9058SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
809999d9058SBarry Smith       req_size[idex] = 0;
810999d9058SBarry Smith       rbuf1_i         = rbuf1[idex];
811d5b485abSSatish Balay       start           = 2*rbuf1_i[0] + 1;
812b24ad042SBarry Smith       ierr            = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
813b24ad042SBarry Smith       ierr            = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr);
814999d9058SBarry Smith       sbuf2_i         = sbuf2[idex];
815d5b485abSSatish Balay       for (j=start; j<end; j++) {
816d5b485abSSatish Balay         id               = rbuf1_i[j] - rstart;
817d5b485abSSatish Balay         ncols            = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
818d5b485abSSatish Balay         sbuf2_i[j]       = ncols;
819999d9058SBarry Smith         req_size[idex] += ncols;
820d5b485abSSatish Balay       }
821999d9058SBarry Smith       req_source[idex] = r_status1[i].MPI_SOURCE;
822d5b485abSSatish Balay       /* form the header */
823999d9058SBarry Smith       sbuf2_i[0]   = req_size[idex];
824d5b485abSSatish Balay       for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
825b24ad042SBarry Smith       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
826d5b485abSSatish Balay     }
827d5b485abSSatish Balay   }
828606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
829606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
830d5b485abSSatish Balay 
831d5b485abSSatish Balay   /*  recv buffer sizes */
832d5b485abSSatish Balay   /* Receive messages*/
833d5b485abSSatish Balay 
834b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr);
835b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr);
836b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr);
837b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr);
838b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr);
839d5b485abSSatish Balay 
840d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
841999d9058SBarry Smith     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
842b24ad042SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr);
843999d9058SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr);
844b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
845b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
846d5b485abSSatish Balay   }
847606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
848606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
849d5b485abSSatish Balay 
850d5b485abSSatish Balay   /* Wait on sends1 and sends2 */
851b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr);
852b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr);
853d5b485abSSatish Balay 
8540c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
8550c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
856606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
857606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
858606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
859606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
860d5b485abSSatish Balay 
861d5b485abSSatish Balay   /* Now allocate buffers for a->j, and send them off */
862b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr);
863d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
864b24ad042SBarry Smith   ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr);
865d5b485abSSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
866d5b485abSSatish Balay 
867b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr);
868d5b485abSSatish Balay   {
869d5b485abSSatish Balay      for (i=0; i<nrqr; i++) {
870d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
871d5b485abSSatish Balay       sbuf_aj_i = sbuf_aj[i];
872d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0] + 1;
873d5b485abSSatish Balay       ct2       = 0;
874d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
875d5b485abSSatish Balay         kmax = rbuf1[i][2*j];
876d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
877e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
878d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
879d5b485abSSatish Balay           ncols  = nzA + nzB;
880d5b485abSSatish Balay           cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
881d5b485abSSatish Balay 
882d5b485abSSatish Balay           /* load the column indices for this row into cols*/
883d5b485abSSatish Balay           cols  = sbuf_aj_i + ct2;
884d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
885d5b485abSSatish Balay             if ((ctmp = bmap[cworkB[l]]) < cstart)  cols[l] = ctmp;
886d5b485abSSatish Balay             else break;
887d5b485abSSatish Balay           }
888d5b485abSSatish Balay           imark = l;
889d5b485abSSatish Balay           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
890d5b485abSSatish Balay           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
891d5b485abSSatish Balay           ct2 += ncols;
892d5b485abSSatish Balay         }
893d5b485abSSatish Balay       }
894b24ad042SBarry Smith       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
895d5b485abSSatish Balay     }
896d5b485abSSatish Balay   }
897b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr);
898b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr);
899d5b485abSSatish Balay 
900d5b485abSSatish Balay   /* Allocate buffers for a->a, and send them off */
90182502324SSatish Balay   ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr);
902d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
90382502324SSatish Balay   ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr);
904a2feddc5SSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
905d5b485abSSatish Balay 
906b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr);
907d5b485abSSatish Balay   {
908d5b485abSSatish Balay     for (i=0; i<nrqr; i++) {
909d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
910d5b485abSSatish Balay       sbuf_aa_i = sbuf_aa[i];
911d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0]+1;
912d5b485abSSatish Balay       ct2       = 0;
913d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
914d5b485abSSatish Balay         kmax = rbuf1_i[2*j];
915d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
916e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
917d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
918d5b485abSSatish Balay           ncols  = nzA + nzB;
919d5b485abSSatish Balay           cworkA = a_j + a_i[row];     cworkB = b_j + b_i[row];
920a2feddc5SSatish Balay           vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2;
921d5b485abSSatish Balay 
922d5b485abSSatish Balay           /* load the column values for this row into vals*/
9235b83ace0SSatish Balay           vals  = sbuf_aa_i+ct2*bs2;
924d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
9253a40ed3dSBarry Smith             if ((bmap[cworkB[l]]) < cstart) {
9263eda8832SBarry Smith               ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9273a40ed3dSBarry Smith             }
928d5b485abSSatish Balay             else break;
929d5b485abSSatish Balay           }
930d5b485abSSatish Balay           imark = l;
9313a40ed3dSBarry Smith           for (l=0; l<nzA; l++) {
9323eda8832SBarry Smith             ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9333a40ed3dSBarry Smith           }
9343a40ed3dSBarry Smith           for (l=imark; l<nzB; l++) {
9353eda8832SBarry Smith             ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9363a40ed3dSBarry Smith           }
937d5b485abSSatish Balay           ct2 += ncols;
938d5b485abSSatish Balay         }
939d5b485abSSatish Balay       }
9403eda8832SBarry Smith       ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
941d5b485abSSatish Balay     }
942d5b485abSSatish Balay   }
943b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr);
944b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr);
945533163c2SBarry Smith   ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr);
946606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
947d5b485abSSatish Balay 
948d5b485abSSatish Balay   /* Form the matrix */
949d5b485abSSatish Balay   /* create col map */
950d5b485abSSatish Balay   {
9515d0c19d7SBarry Smith     const PetscInt *icol_i;
952233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
953233c2ff6SSatish Balay     /* Create row map*/
954b0a32e0cSBarry Smith     ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&colmaps);CHKERRQ(ierr);
955ff0794f7SSatish Balay     for (i=0; i<ismax; i++) {
956ff0794f7SSatish Balay       ierr = PetscTableCreate(ncol[i]+1,&colmaps[i]);CHKERRQ(ierr);
957233c2ff6SSatish Balay     }
958233c2ff6SSatish Balay #else
95923bdbc58SBarry Smith     ierr    = PetscMalloc(ismax*sizeof(PetscInt*),&cmap);CHKERRQ(ierr);
96023bdbc58SBarry Smith     ierr    = PetscMalloc(ismax*c->Nbs*sizeof(PetscInt),&cmap[0]);CHKERRQ(ierr);
961a2feddc5SSatish Balay     for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; }
962233c2ff6SSatish Balay #endif
963d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
964d5b485abSSatish Balay       jmax   = ncol[i];
965d5b485abSSatish Balay       icol_i = icol[i];
966233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
967233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
968233c2ff6SSatish Balay       for (j=0; j<jmax; j++) {
969001aedefSBarry Smith 	ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr);
970233c2ff6SSatish Balay       }
971233c2ff6SSatish Balay #else
972d5b485abSSatish Balay       cmap_i = cmap[i];
973d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
974d5b485abSSatish Balay         cmap_i[icol_i[j]] = j+1;
975d5b485abSSatish Balay       }
976233c2ff6SSatish Balay #endif
977d5b485abSSatish Balay     }
978d5b485abSSatish Balay   }
979d5b485abSSatish Balay 
980d5b485abSSatish Balay   /* Create lens which is required for MatCreate... */
981d5b485abSSatish Balay   for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
982052f0c41SBarry Smith   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr);
983b24ad042SBarry Smith   lens[0] = (PetscInt*)(lens + ismax);
984b24ad042SBarry Smith   ierr    = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr);
985d5b485abSSatish Balay   for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
986d5b485abSSatish Balay 
987d5b485abSSatish Balay   /* Update lens from local data */
988d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
989d5b485abSSatish Balay     jmax   = nrow[i];
990233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
991233c2ff6SSatish Balay     lcol1_gcol1 = colmaps[i];
992233c2ff6SSatish Balay #else
993d5b485abSSatish Balay     cmap_i = cmap[i];
994233c2ff6SSatish Balay #endif
995d5b485abSSatish Balay     irow_i = irow[i];
996d5b485abSSatish Balay     lens_i = lens[i];
997d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
998d5b485abSSatish Balay       row  = irow_i[j];
999233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1000899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1001233c2ff6SSatish Balay #else
1002d5b485abSSatish Balay       proc = rtable[row];
1003233c2ff6SSatish Balay #endif
1004d5b485abSSatish Balay       if (proc == rank) {
100536f4e84dSSatish Balay         /* Get indices from matA and then from matB */
100636f4e84dSSatish Balay         row    = row - rstart;
100736f4e84dSSatish Balay         nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
100836f4e84dSSatish Balay         cworkA =  a_j + a_i[row]; cworkB = b_j + b_i[row];
1009233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1010233c2ff6SSatish Balay        for (k=0; k<nzA; k++) {
1011233c2ff6SSatish Balay 	 ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1012233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1013233c2ff6SSatish Balay         }
1014233c2ff6SSatish Balay         for (k=0; k<nzB; k++) {
1015233c2ff6SSatish Balay 	  ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1016233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1017233c2ff6SSatish Balay         }
1018233c2ff6SSatish Balay #else
1019ca161407SBarry Smith         for (k=0; k<nzA; k++) {
102036f4e84dSSatish Balay           if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; }
1021ca161407SBarry Smith         }
1022ca161407SBarry Smith         for (k=0; k<nzB; k++) {
102336f4e84dSSatish Balay           if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; }
1024d5b485abSSatish Balay         }
1025233c2ff6SSatish Balay #endif
1026d5b485abSSatish Balay       }
1027a2feddc5SSatish Balay     }
1028ca161407SBarry Smith   }
1029233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1030d5b485abSSatish Balay   /* Create row map*/
1031b0a32e0cSBarry Smith   ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rowmaps);CHKERRQ(ierr);
1032ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
1033ff0794f7SSatish Balay     ierr = PetscTableCreate(nrow[i]+1,&rowmaps[i]);CHKERRQ(ierr);
1034233c2ff6SSatish Balay   }
1035233c2ff6SSatish Balay #else
1036233c2ff6SSatish Balay   /* Create row map*/
1037052f0c41SBarry Smith   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr);
1038b24ad042SBarry Smith   rmap[0] = (PetscInt*)(rmap + ismax);
1039b24ad042SBarry Smith   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr);
1040233c2ff6SSatish Balay   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;}
1041233c2ff6SSatish Balay #endif
1042d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1043d5b485abSSatish Balay     irow_i = irow[i];
1044d5b485abSSatish Balay     jmax   = nrow[i];
1045233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1046233c2ff6SSatish Balay     lrow1_grow1 = rowmaps[i];
1047233c2ff6SSatish Balay     for (j=0; j<jmax; j++) {
1048233c2ff6SSatish Balay       ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr);
1049233c2ff6SSatish Balay     }
1050233c2ff6SSatish Balay #else
1051233c2ff6SSatish Balay     rmap_i = rmap[i];
1052d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1053d5b485abSSatish Balay       rmap_i[irow_i[j]] = j;
1054d5b485abSSatish Balay     }
1055233c2ff6SSatish Balay #endif
1056d5b485abSSatish Balay   }
1057d5b485abSSatish Balay 
1058d5b485abSSatish Balay   /* Update lens from offproc data */
1059d5b485abSSatish Balay   {
1060b24ad042SBarry Smith     PetscInt    *rbuf2_i,*rbuf3_i,*sbuf1_i;
1061b24ad042SBarry Smith     PetscMPIInt ii;
1062d5b485abSSatish Balay 
1063d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1064b24ad042SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr);
1065b24ad042SBarry Smith       idex   = pa[ii];
1066999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1067d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1068d5b485abSSatish Balay       ct1     = 2*jmax+1;
1069d5b485abSSatish Balay       ct2     = 0;
1070b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1071b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
1072d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1073d5b485abSSatish Balay         is_no   = sbuf1_i[2*j-1];
1074d5b485abSSatish Balay         max1    = sbuf1_i[2*j];
1075d5b485abSSatish Balay         lens_i  = lens[is_no];
1076233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1077233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1078233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1079233c2ff6SSatish Balay #else
1080d5b485abSSatish Balay         cmap_i  = cmap[is_no];
1081d5b485abSSatish Balay 	rmap_i  = rmap[is_no];
1082233c2ff6SSatish Balay #endif
1083d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1084233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1085233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1086233c2ff6SSatish Balay           row--;
1087e32f2f54SBarry Smith           if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1088233c2ff6SSatish Balay #else
1089d5b485abSSatish Balay           row  = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
1090233c2ff6SSatish Balay #endif
1091d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1092d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1093233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1094233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
1095233c2ff6SSatish Balay 	    if (tt) {
1096233c2ff6SSatish Balay               lens_i[row]++;
1097233c2ff6SSatish Balay             }
1098233c2ff6SSatish Balay #else
1099d5b485abSSatish Balay             if (cmap_i[rbuf3_i[ct2]]) {
1100d5b485abSSatish Balay               lens_i[row]++;
1101d5b485abSSatish Balay             }
1102233c2ff6SSatish Balay #endif
1103d5b485abSSatish Balay           }
1104d5b485abSSatish Balay         }
1105d5b485abSSatish Balay       }
1106d5b485abSSatish Balay     }
1107d5b485abSSatish Balay   }
1108606d414cSSatish Balay   ierr = PetscFree(r_status3);CHKERRQ(ierr);
1109606d414cSSatish Balay   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
11100c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
1111606d414cSSatish Balay   ierr = PetscFree(s_status3);CHKERRQ(ierr);
1112606d414cSSatish Balay   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1113d5b485abSSatish Balay 
1114d5b485abSSatish Balay   /* Create the submatrices */
1115d5b485abSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1116d5b485abSSatish Balay     /*
1117d5b485abSSatish Balay         Assumes new rows are same length as the old rows, hence bug!
1118d5b485abSSatish Balay     */
1119d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1120df36731bSBarry Smith       mat = (Mat_SeqBAIJ *)(submats[i]->data);
1121e7e72b3dSBarry Smith       if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1122b24ad042SBarry Smith       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
1123e7e72b3dSBarry Smith       if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1124d5b485abSSatish Balay       /* Initial matrix as if empty */
1125b24ad042SBarry Smith       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr);
1126d5f3da31SBarry Smith       submats[i]->factortype = C->factortype;
1127d5b485abSSatish Balay     }
1128ca161407SBarry Smith   } else {
1129d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1130f69a0ea3SMatthew Knepley       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
1131f69a0ea3SMatthew Knepley       ierr = MatSetSizes(submats[i],nrow[i]*bs,ncol[i]*bs,nrow[i]*bs,ncol[i]*bs);CHKERRQ(ierr);
11327adad957SLisandro Dalcin       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
1133d0f46423SBarry Smith       ierr = MatSeqBAIJSetPreallocation(submats[i],C->rmap->bs,0,lens[i]);CHKERRQ(ierr);
1134a6d3ce06SHong Zhang       ierr = MatSeqSBAIJSetPreallocation(submats[i],C->rmap->bs,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by MatGetSubMatrices_MPISBAIJ()*/
1135d5b485abSSatish Balay     }
1136d5b485abSSatish Balay   }
1137d5b485abSSatish Balay 
1138d5b485abSSatish Balay   /* Assemble the matrices */
1139d5b485abSSatish Balay   /* First assemble the local rows */
1140d5b485abSSatish Balay   {
1141b24ad042SBarry Smith     PetscInt       ilen_row,*imat_ilen,*imat_j,*imat_i;
11423eda8832SBarry Smith     MatScalar *imat_a;
1143d5b485abSSatish Balay 
1144d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1145df36731bSBarry Smith       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1146d5b485abSSatish Balay       imat_ilen = mat->ilen;
1147d5b485abSSatish Balay       imat_j    = mat->j;
1148d5b485abSSatish Balay       imat_i    = mat->i;
1149d5b485abSSatish Balay       imat_a    = mat->a;
1150233c2ff6SSatish Balay 
1151233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1152233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
1153233c2ff6SSatish Balay       lrow1_grow1 = rowmaps[i];
1154233c2ff6SSatish Balay #else
1155d5b485abSSatish Balay       cmap_i    = cmap[i];
1156d5b485abSSatish Balay       rmap_i    = rmap[i];
1157233c2ff6SSatish Balay #endif
1158d5b485abSSatish Balay       irow_i    = irow[i];
1159d5b485abSSatish Balay       jmax      = nrow[i];
1160d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1161d5b485abSSatish Balay         row      = irow_i[j];
1162233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1163899cda47SBarry Smith 	ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1164233c2ff6SSatish Balay #else
1165d5b485abSSatish Balay 	proc = rtable[row];
1166233c2ff6SSatish Balay #endif
1167d5b485abSSatish Balay         if (proc == rank) {
116836f4e84dSSatish Balay           row      = row - rstart;
116936f4e84dSSatish Balay           nzA      = a_i[row+1] - a_i[row];
117036f4e84dSSatish Balay           nzB      = b_i[row+1] - b_i[row];
117136f4e84dSSatish Balay           cworkA   = a_j + a_i[row];
117236f4e84dSSatish Balay           cworkB   = b_j + b_i[row];
117336f4e84dSSatish Balay           vworkA   = a_a + a_i[row]*bs2;
117436f4e84dSSatish Balay           vworkB   = b_a + b_i[row]*bs2;
1175233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1176233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr);
1177233c2ff6SSatish Balay           row--;
1178e32f2f54SBarry Smith           if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1179233c2ff6SSatish Balay #else
118036f4e84dSSatish Balay           row      = rmap_i[row + rstart];
1181233c2ff6SSatish Balay #endif
1182df36731bSBarry Smith           mat_i    = imat_i[row];
118336f4e84dSSatish Balay           mat_a    = imat_a + mat_i*bs2;
1184d5b485abSSatish Balay           mat_j    = imat_j + mat_i;
118536f4e84dSSatish Balay           ilen_row = imat_ilen[row];
118636f4e84dSSatish Balay 
118736f4e84dSSatish Balay           /* load the column indices for this row into cols*/
118836f4e84dSSatish Balay           for (l=0; l<nzB; l++) {
118936f4e84dSSatish Balay 	    if ((ctmp = bmap[cworkB[l]]) < cstart) {
1190233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1191233c2ff6SSatish Balay 	      ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr);
1192233c2ff6SSatish Balay 	      if (tcol) {
1193233c2ff6SSatish Balay #else
119436f4e84dSSatish Balay               if ((tcol = cmap_i[ctmp])) {
1195233c2ff6SSatish Balay #endif
1196df36731bSBarry Smith                 *mat_j++ = tcol - 1;
11973eda8832SBarry Smith                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1198549d3d68SSatish Balay                 mat_a   += bs2;
1199d5b485abSSatish Balay                 ilen_row++;
1200d5b485abSSatish Balay               }
1201ca161407SBarry Smith             } else break;
120236f4e84dSSatish Balay           }
120336f4e84dSSatish Balay           imark = l;
120436f4e84dSSatish Balay           for (l=0; l<nzA; l++) {
1205233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1206233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1207233c2ff6SSatish Balay 	    if (tcol) {
1208233c2ff6SSatish Balay #else
120936f4e84dSSatish Balay 	    if ((tcol = cmap_i[cstart + cworkA[l]])) {
1210233c2ff6SSatish Balay #endif
121136f4e84dSSatish Balay               *mat_j++ = tcol - 1;
12123eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1213549d3d68SSatish Balay               mat_a   += bs2;
121436f4e84dSSatish Balay               ilen_row++;
121536f4e84dSSatish Balay             }
121636f4e84dSSatish Balay           }
121736f4e84dSSatish Balay           for (l=imark; l<nzB; l++) {
1218233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1219233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1220233c2ff6SSatish Balay 	    if (tcol) {
1221233c2ff6SSatish Balay #else
122236f4e84dSSatish Balay             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1223233c2ff6SSatish Balay #endif
122436f4e84dSSatish Balay               *mat_j++ = tcol - 1;
12253eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1226549d3d68SSatish Balay               mat_a   += bs2;
122736f4e84dSSatish Balay               ilen_row++;
122836f4e84dSSatish Balay             }
122936f4e84dSSatish Balay           }
1230d5b485abSSatish Balay           imat_ilen[row] = ilen_row;
1231d5b485abSSatish Balay         }
1232d5b485abSSatish Balay       }
123336f4e84dSSatish Balay 
1234d5b485abSSatish Balay     }
1235d5b485abSSatish Balay   }
1236d5b485abSSatish Balay 
1237d5b485abSSatish Balay   /*   Now assemble the off proc rows*/
1238d5b485abSSatish Balay   {
1239b24ad042SBarry Smith     PetscInt    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1240b24ad042SBarry Smith     PetscInt    *imat_j,*imat_i;
12413eda8832SBarry Smith     MatScalar   *imat_a,*rbuf4_i;
1242b24ad042SBarry Smith     PetscMPIInt ii;
1243d5b485abSSatish Balay 
1244d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1245b24ad042SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr);
1246b24ad042SBarry Smith       idex   = pa[ii];
1247999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1248d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1249d5b485abSSatish Balay       ct1     = 2*jmax + 1;
1250d5b485abSSatish Balay       ct2     = 0;
1251b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1252b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
1253b24ad042SBarry Smith       rbuf4_i = rbuf4[ii];
1254d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1255d5b485abSSatish Balay         is_no     = sbuf1_i[2*j-1];
1256233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1257233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1258233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1259233c2ff6SSatish Balay #else
1260d5b485abSSatish Balay         rmap_i    = rmap[is_no];
1261d5b485abSSatish Balay         cmap_i    = cmap[is_no];
1262233c2ff6SSatish Balay #endif
1263df36731bSBarry Smith         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1264d5b485abSSatish Balay         imat_ilen = mat->ilen;
1265d5b485abSSatish Balay         imat_j    = mat->j;
1266d5b485abSSatish Balay         imat_i    = mat->i;
1267d5b485abSSatish Balay         imat_a    = mat->a;
1268d5b485abSSatish Balay         max1      = sbuf1_i[2*j];
1269d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1270d5b485abSSatish Balay           row   = sbuf1_i[ct1];
1271233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1272233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr);
1273233c2ff6SSatish Balay           row--;
1274e32f2f54SBarry Smith           if(row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1275233c2ff6SSatish Balay #else
1276d5b485abSSatish Balay           row   = rmap_i[row];
1277233c2ff6SSatish Balay #endif
1278d5b485abSSatish Balay           ilen  = imat_ilen[row];
1279df36731bSBarry Smith           mat_i = imat_i[row];
128036f4e84dSSatish Balay           mat_a = imat_a + mat_i*bs2;
1281d5b485abSSatish Balay           mat_j = imat_j + mat_i;
1282d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1283d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1284233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1285233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1286233c2ff6SSatish Balay 	    if (tcol) {
1287233c2ff6SSatish Balay #else
1288d5b485abSSatish Balay 	    if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1289233c2ff6SSatish Balay #endif
1290df36731bSBarry Smith               *mat_j++    = tcol - 1;
129136f4e84dSSatish Balay               /* *mat_a++= rbuf4_i[ct2]; */
12923eda8832SBarry Smith               ierr        = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1293549d3d68SSatish Balay               mat_a      += bs2;
1294d5b485abSSatish Balay               ilen++;
1295d5b485abSSatish Balay             }
1296d5b485abSSatish Balay           }
1297d5b485abSSatish Balay           imat_ilen[row] = ilen;
1298d5b485abSSatish Balay         }
1299d5b485abSSatish Balay       }
1300d5b485abSSatish Balay     }
1301d5b485abSSatish Balay   }
1302606d414cSSatish Balay   ierr = PetscFree(r_status4);CHKERRQ(ierr);
1303606d414cSSatish Balay   ierr = PetscFree(r_waits4);CHKERRQ(ierr);
13040c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
1305606d414cSSatish Balay   ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1306606d414cSSatish Balay   ierr = PetscFree(s_status4);CHKERRQ(ierr);
1307d5b485abSSatish Balay 
1308d5b485abSSatish Balay   /* Restore the indices */
1309d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1310d5b485abSSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1311d5b485abSSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1312d5b485abSSatish Balay   }
1313d5b485abSSatish Balay 
1314d5b485abSSatish Balay   /* Destroy allocated memory */
131523bdbc58SBarry Smith #if defined(PETSC_USE_CTABLE)
131623bdbc58SBarry Smith   ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr);
131723bdbc58SBarry Smith #else
131823bdbc58SBarry Smith   ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr);
131923bdbc58SBarry Smith #endif
132023bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
1321606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
1322d5b485abSSatish Balay 
132323bdbc58SBarry Smith   ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
1324606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1325606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1326d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
1327606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1328d5b485abSSatish Balay   }
1329d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1330606d414cSSatish Balay     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1331606d414cSSatish Balay     ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);
1332d5b485abSSatish Balay   }
133323bdbc58SBarry Smith   ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr);
1334606d414cSSatish Balay   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1335606d414cSSatish Balay   ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1336606d414cSSatish Balay   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1337606d414cSSatish Balay   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
1338606d414cSSatish Balay   ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1339606d414cSSatish Balay   ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
1340d5b485abSSatish Balay 
1341233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1342ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
13439c666560SBarry Smith     ierr = PetscTableDestroy(rowmaps[i]);CHKERRQ(ierr);
13449c666560SBarry Smith     ierr = PetscTableDestroy(colmaps[i]);CHKERRQ(ierr);
1345233c2ff6SSatish Balay   }
1346233c2ff6SSatish Balay   ierr = PetscFree(colmaps);CHKERRQ(ierr);
1347233c2ff6SSatish Balay   ierr = PetscFree(rowmaps);CHKERRQ(ierr);
1348233c2ff6SSatish Balay #else
1349606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
135023bdbc58SBarry Smith   ierr = PetscFree(cmap[0]);CHKERRQ(ierr);
1351233c2ff6SSatish Balay   ierr = PetscFree(cmap);CHKERRQ(ierr);
1352233c2ff6SSatish Balay #endif
1353606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1354d5b485abSSatish Balay 
1355d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
135636f4e84dSSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
135736f4e84dSSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1358d5b485abSSatish Balay   }
13593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1360d5b485abSSatish Balay }
1361ca161407SBarry Smith 
1362