xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 4da72fa9edbadbed037b02b2107ab8f0d233e822)
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 */
2505d8c843SHong Zhang   ierr = ISCompressIndicesGeneral(N,C->rmap->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   }
306bf464f9SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);}
3105d8c843SHong Zhang   ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr);
326bf464f9SBarry 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"
62db41cccfSHong 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) {
2336bf464f9SBarry Smith     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 
5384a2ae208SSatish Balay #undef __FUNCT__
5394a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ"
540b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
541d5b485abSSatish Balay {
54236f4e84dSSatish Balay   IS             *isrow_new,*iscol_new;
543cf4f527aSSatish Balay   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
5446849ba73SBarry Smith   PetscErrorCode ierr;
545*4da72fa9SHong Zhang   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs;
546*4da72fa9SHong Zhang   PetscBool      colflag,*allcolumns,*allrows;
547a2feddc5SSatish Balay 
5483a40ed3dSBarry Smith   PetscFunctionBegin;
5493a40ed3dSBarry Smith   /* The compression and expansion should be avoided. Does'nt point
550a2feddc5SSatish Balay      out errors might change the indices hence buggey */
551a2feddc5SSatish Balay 
55223bdbc58SBarry Smith   ierr = PetscMalloc2(ismax+1,IS,&isrow_new,ismax+1,IS,&iscol_new);CHKERRQ(ierr);
55305d8c843SHong Zhang   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_new);CHKERRQ(ierr);
55405d8c843SHong Zhang   ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_new);CHKERRQ(ierr);
555cf4f527aSSatish Balay 
556307b7a18SHong Zhang   /* Check for special case: each processor gets entire matrix columns */
557*4da72fa9SHong Zhang   ierr = PetscMalloc2(ismax+1,PetscBool,&allcolumns,ismax+1,PetscBool,&allrows);CHKERRQ(ierr);
558307b7a18SHong Zhang   for (i=0; i<ismax; i++) {
559307b7a18SHong Zhang     ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr);
560307b7a18SHong Zhang     ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr);
561307b7a18SHong Zhang     if (colflag && ncol == C->cmap->N){
562307b7a18SHong Zhang       allcolumns[i] = PETSC_TRUE;
563307b7a18SHong Zhang     } else {
564307b7a18SHong Zhang       allcolumns[i] = PETSC_FALSE;
565307b7a18SHong Zhang     }
566*4da72fa9SHong Zhang 
567*4da72fa9SHong Zhang     ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr);
568*4da72fa9SHong Zhang     ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr);
569*4da72fa9SHong Zhang     if (colflag && nrow == C->rmap->N){
570*4da72fa9SHong Zhang       allrows[i] = PETSC_TRUE;
571*4da72fa9SHong Zhang     } else {
572*4da72fa9SHong Zhang       allrows[i] = PETSC_FALSE;
573*4da72fa9SHong Zhang     }
574307b7a18SHong Zhang   }
575307b7a18SHong Zhang 
576cf4f527aSSatish Balay   /* Allocate memory to hold all the submatrices */
577cf4f527aSSatish Balay   if (scall != MAT_REUSE_MATRIX) {
57882502324SSatish Balay     ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr);
579cf4f527aSSatish Balay   }
580cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
581b24ad042SBarry Smith   nmax          = 20*1000000 / (c->Nbs * sizeof(PetscInt));
582329f5518SBarry Smith   if (!nmax) nmax = 1;
583cf4f527aSSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
584cf4f527aSSatish Balay 
585653e4784SBarry Smith   /* Make sure every processor loops through the nstages */
5867adad957SLisandro Dalcin   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);CHKERRQ(ierr);
587cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
588cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
589cf4f527aSSatish Balay     else if (pos == ismax) max_no = 0;
590cf4f527aSSatish Balay     else                   max_no = ismax-pos;
591*4da72fa9SHong Zhang     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr);
592cf4f527aSSatish Balay     pos += max_no;
593cf4f527aSSatish Balay   }
59436f4e84dSSatish Balay 
59536f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
5966bf464f9SBarry Smith     ierr = ISDestroy(&isrow_new[i]);CHKERRQ(ierr);
5976bf464f9SBarry Smith     ierr = ISDestroy(&iscol_new[i]);CHKERRQ(ierr);
59836f4e84dSSatish Balay   }
59923bdbc58SBarry Smith   ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr);
600*4da72fa9SHong Zhang   ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr);
6013a40ed3dSBarry Smith   PetscFunctionReturn(0);
602a2feddc5SSatish Balay }
603a2feddc5SSatish Balay 
604233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
6054a2ae208SSatish Balay #undef __FUNCT__
6064a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc"
607e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
608233c2ff6SSatish Balay {
609e005ede5SBarry Smith   PetscInt    nGlobalNd = proc_gnode[size];
610b9f7ace7SBarry Smith   PetscMPIInt fproc = PetscMPIIntCast( (PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)));
611233c2ff6SSatish Balay 
612233c2ff6SSatish Balay   PetscFunctionBegin;
61323ce1328SBarry Smith   if (fproc > size) fproc = size;
614e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
615e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
616233c2ff6SSatish Balay     else                         fproc++;
617233c2ff6SSatish Balay   }
618e005ede5SBarry Smith   *rank = fproc;
619233c2ff6SSatish Balay   PetscFunctionReturn(0);
620233c2ff6SSatish Balay }
621233c2ff6SSatish Balay #endif
622233c2ff6SSatish Balay 
623a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
624b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
6254a2ae208SSatish Balay #undef __FUNCT__
6264a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local"
627*4da72fa9SHong Zhang PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats)
628a2feddc5SSatish Balay {
629df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
630cf4f527aSSatish Balay   Mat            A = c->A;
631df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
6325d0c19d7SBarry Smith   const PetscInt **irow,**icol,*irow_i;
6335d0c19d7SBarry Smith   PetscInt       *nrow,*ncol,*w3,*w4,start;
6346849ba73SBarry Smith   PetscErrorCode ierr;
6359405f653SBarry Smith   PetscMPIInt    size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc;
6369405f653SBarry Smith   PetscInt       **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row;
637b24ad042SBarry Smith   PetscInt       nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol;
638b24ad042SBarry Smith   PetscInt       **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
6395d0c19d7SBarry Smith   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
640052f0c41SBarry Smith   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
641d0f46423SBarry Smith   PetscInt       bs=C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
642899cda47SBarry Smith   PetscInt       cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
643899cda47SBarry Smith   PetscInt       *bmap = c->garray,ctmp,rstart=c->rstartbs;
6447a868f3eSHong Zhang   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3;
6457a868f3eSHong Zhang   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3;
646d5b485abSSatish Balay   MPI_Comm       comm;
647ace3abfcSBarry Smith   PetscBool      flag;
648b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1;
6497a868f3eSHong Zhang   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
6507a868f3eSHong Zhang   /* variables below are used for the matrix numerical values - case of !ijonly */
6517a868f3eSHong Zhang   MPI_Request    *r_waits4,*s_waits4;
6527a868f3eSHong Zhang   MPI_Status     *r_status4,*s_status4;
653f3e30f2cSMatthew G Knepley   MatScalar      **rbuf4,**sbuf_aa,*vals,*mat_a = PETSC_NULL,*sbuf_aa_i,*vworkA = PETSC_NULL,*vworkB = PETSC_NULL;
6547a868f3eSHong Zhang   MatScalar      *a_a=a->a,*b_a=b->a;
655c7dd2924SSatish Balay 
65680d608c0SSatish Balay #if defined (PETSC_USE_CTABLE)
657b24ad042SBarry Smith   PetscInt       tt;
6589066e3d5SHong Zhang   PetscTable     *rmap,*cmap,rmap_i,cmap_i=PETSC_NULL;
659233c2ff6SSatish Balay #else
6609066e3d5SHong Zhang   PetscInt       **cmap,*cmap_i=PETSC_NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
661233c2ff6SSatish Balay #endif
662d5b485abSSatish Balay 
6633a40ed3dSBarry Smith   PetscFunctionBegin;
6647adad957SLisandro Dalcin   comm   = ((PetscObject)C)->comm;
6657adad957SLisandro Dalcin   tag0   = ((PetscObject)C)->tag;
666d5b485abSSatish Balay   size   = c->size;
667d5b485abSSatish Balay   rank   = c->rank;
668d5b485abSSatish Balay 
66934e6ae18SSatish Balay   /* Get some new tags to keep the communication clean */
67034e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
67134e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
67234e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
67334e6ae18SSatish Balay 
674052f0c41SBarry Smith #if defined(PETSC_USE_CTABLE)
67523bdbc58SBarry Smith   ierr = PetscMalloc4(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol);CHKERRQ(ierr);
676052f0c41SBarry Smith #else
67723bdbc58SBarry Smith   ierr = PetscMalloc5(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol,Mbs+1,PetscInt,&rtable);CHKERRQ(ierr);
678d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
679d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
68082a7e548SBarry Smith     jmax = C->rmap->range[i+1]/bs;
681d5b485abSSatish Balay     for (; j<jmax; j++) {
682d5b485abSSatish Balay       rtable[j] = i;
683d5b485abSSatish Balay     }
684d5b485abSSatish Balay   }
685233c2ff6SSatish Balay #endif
686233c2ff6SSatish Balay 
687233c2ff6SSatish Balay   for (i=0; i<ismax; i++) {
688*4da72fa9SHong Zhang     if (allrows[i]){
689*4da72fa9SHong Zhang       irow[i] = PETSC_NULL;
690*4da72fa9SHong Zhang       nrow[i] = C->rmap->N/bs;
691*4da72fa9SHong Zhang     } else {
692233c2ff6SSatish Balay       ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
693233c2ff6SSatish Balay       ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
694*4da72fa9SHong Zhang     }
695*4da72fa9SHong Zhang 
696307b7a18SHong Zhang     if (allcolumns[i]){
697307b7a18SHong Zhang       icol[i] = PETSC_NULL;
698307b7a18SHong Zhang       ncol[i] = C->cmap->N/bs;
699307b7a18SHong Zhang     } else {
700307b7a18SHong Zhang       ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
701233c2ff6SSatish Balay       ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
702233c2ff6SSatish Balay     }
703307b7a18SHong Zhang   }
704d5b485abSSatish Balay 
705d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg,and buffer space
706d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
70723bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr);
70823bdbc58SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
70923bdbc58SBarry Smith   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
71023bdbc58SBarry Smith   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
711d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
712b24ad042SBarry Smith     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
713d5b485abSSatish Balay     jmax   = nrow[i];
714d5b485abSSatish Balay     irow_i = irow[i];
715d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
716*4da72fa9SHong Zhang       if (allrows[i]) {
717*4da72fa9SHong Zhang         row = j;
718*4da72fa9SHong Zhang       } else {
719d5b485abSSatish Balay         row  = irow_i[j];
720*4da72fa9SHong Zhang       }
721233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
722899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
723233c2ff6SSatish Balay #else
724d5b485abSSatish Balay       proc = rtable[row];
725233c2ff6SSatish Balay #endif
726d5b485abSSatish Balay       w4[proc]++;
727d5b485abSSatish Balay     }
728d5b485abSSatish Balay     for (j=0; j<size; j++) {
729d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
730d5b485abSSatish Balay     }
731d5b485abSSatish Balay   }
732d5b485abSSatish Balay 
733d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
734e757655aSSatish Balay   msz      = 0;              /* total mesg length for all proc */
735d5b485abSSatish Balay   w1[rank] = 0;              /* no mesg sent to intself */
736d5b485abSSatish Balay   w3[rank] = 0;
737d5b485abSSatish Balay   for (i=0; i<size; i++) {
738d5b485abSSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
739d5b485abSSatish Balay   }
740b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/
741d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
742d5b485abSSatish Balay     if (w1[i]) { pa[j] = i; j++; }
743d5b485abSSatish Balay   }
744d5b485abSSatish Balay 
745d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
746d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
747d5b485abSSatish Balay     j     = pa[i];
748d5b485abSSatish Balay     w1[j] += w2[j] + 2* w3[j];
749d5b485abSSatish Balay     msz   += w1[j];
750d5b485abSSatish Balay   }
751d5b485abSSatish Balay 
752c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
753a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
754c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
755d5b485abSSatish Balay 
756c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
757c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
758c7dd2924SSatish Balay 
759c7dd2924SSatish Balay   ierr = PetscFree(onodes1);CHKERRQ(ierr);
760c7dd2924SSatish Balay   ierr = PetscFree(olengths1);CHKERRQ(ierr);
761d5b485abSSatish Balay 
762d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
76323bdbc58SBarry Smith   ierr = PetscMalloc4(size,PetscInt*,&sbuf1,size,PetscInt*,&ptr,2*msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr);
76423bdbc58SBarry Smith   ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
76523bdbc58SBarry Smith   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
766d5b485abSSatish Balay   {
767b24ad042SBarry Smith     PetscInt *iptr = tmp,ict = 0;
768d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
769d5b485abSSatish Balay       j         = pa[i];
770d5b485abSSatish Balay       iptr     += ict;
771d5b485abSSatish Balay       sbuf1[j]  = iptr;
772d5b485abSSatish Balay       ict       = w1[j];
773d5b485abSSatish Balay     }
774d5b485abSSatish Balay   }
775d5b485abSSatish Balay 
776d5b485abSSatish Balay   /* Form the outgoing messages */
777d5b485abSSatish Balay   /* Initialise the header space */
778d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
779d5b485abSSatish Balay     j           = pa[i];
780d5b485abSSatish Balay     sbuf1[j][0] = 0;
781b24ad042SBarry Smith     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
782d5b485abSSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
783d5b485abSSatish Balay   }
784d5b485abSSatish Balay 
785d5b485abSSatish Balay   /* Parse the isrow and copy data into outbuf */
786d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
787b24ad042SBarry Smith     ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
788d5b485abSSatish Balay     irow_i = irow[i];
789d5b485abSSatish Balay     jmax   = nrow[i];
790d5b485abSSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
791*4da72fa9SHong Zhang       if (allrows[i]){
792*4da72fa9SHong Zhang         row = j;
793*4da72fa9SHong Zhang       } else {
794d5b485abSSatish Balay         row  = irow_i[j];
795*4da72fa9SHong Zhang       }
796233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
797899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
798233c2ff6SSatish Balay #else
799d5b485abSSatish Balay       proc = rtable[row];
800233c2ff6SSatish Balay #endif
801d5b485abSSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
802d5b485abSSatish Balay         ctr[proc]++;
803d5b485abSSatish Balay         *ptr[proc] = row;
804d5b485abSSatish Balay         ptr[proc]++;
805d5b485abSSatish Balay       }
806d5b485abSSatish Balay     }
807d5b485abSSatish Balay     /* Update the headers for the current IS */
808d5b485abSSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
80906ef35abSSatish Balay       if ((ctr_j = ctr[j])) {
810d5b485abSSatish Balay         sbuf1_j        = sbuf1[j];
811d5b485abSSatish Balay         k              = ++sbuf1_j[0];
812d5b485abSSatish Balay         sbuf1_j[2*k]   = ctr_j;
813d5b485abSSatish Balay         sbuf1_j[2*k-1] = i;
814d5b485abSSatish Balay       }
815d5b485abSSatish Balay     }
816d5b485abSSatish Balay   }
817d5b485abSSatish Balay 
818d5b485abSSatish Balay   /*  Now  post the sends */
819b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
820d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
821d5b485abSSatish Balay     j = pa[i];
822b24ad042SBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
823d5b485abSSatish Balay   }
824d5b485abSSatish Balay 
825d5b485abSSatish Balay   /* Post Recieves to capture the buffer size */
826b0a32e0cSBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr);
827b24ad042SBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr);
828d5b485abSSatish Balay   rbuf2[0] = tmp + msz;
829d5b485abSSatish Balay   for (i=1; i<nrqs; ++i) {
830d5b485abSSatish Balay     j        = pa[i];
831d5b485abSSatish Balay     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
832d5b485abSSatish Balay   }
833d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
834d5b485abSSatish Balay     j    = pa[i];
835b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
836d5b485abSSatish Balay   }
837d5b485abSSatish Balay 
838d5b485abSSatish Balay   /* Send to other procs the buf size they should allocate */
839d5b485abSSatish Balay 
840d5b485abSSatish Balay   /* Receive messages*/
841b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
842b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr);
84323bdbc58SBarry Smith   ierr = PetscMalloc3(nrqr+1,PetscInt*,&sbuf2,nrqr,PetscInt,&req_size,nrqr,PetscInt,&req_source);CHKERRQ(ierr);
844d5b485abSSatish Balay   {
845df36731bSBarry Smith     Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
846b24ad042SBarry Smith     PetscInt    *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
847d5b485abSSatish Balay 
848d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
849999d9058SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
850999d9058SBarry Smith       req_size[idex] = 0;
851999d9058SBarry Smith       rbuf1_i         = rbuf1[idex];
852d5b485abSSatish Balay       start           = 2*rbuf1_i[0] + 1;
853b24ad042SBarry Smith       ierr            = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
854b24ad042SBarry Smith       ierr            = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr);
855999d9058SBarry Smith       sbuf2_i         = sbuf2[idex];
856d5b485abSSatish Balay       for (j=start; j<end; j++) {
857d5b485abSSatish Balay         id               = rbuf1_i[j] - rstart;
858d5b485abSSatish Balay         ncols            = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
859d5b485abSSatish Balay         sbuf2_i[j]       = ncols;
860999d9058SBarry Smith         req_size[idex] += ncols;
861d5b485abSSatish Balay       }
862999d9058SBarry Smith       req_source[idex] = r_status1[i].MPI_SOURCE;
863d5b485abSSatish Balay       /* form the header */
864999d9058SBarry Smith       sbuf2_i[0]   = req_size[idex];
865d5b485abSSatish Balay       for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
866b24ad042SBarry Smith       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
867d5b485abSSatish Balay     }
868d5b485abSSatish Balay   }
869606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
870606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
871d5b485abSSatish Balay 
872d5b485abSSatish Balay   /*  recv buffer sizes */
873d5b485abSSatish Balay   /* Receive messages*/
874b24ad042SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr);
875b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr);
876b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr);
8777a868f3eSHong Zhang   if (!ijonly){
8787a868f3eSHong Zhang     ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr);
8797a868f3eSHong Zhang     ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr);
8807a868f3eSHong Zhang   }
881d5b485abSSatish Balay 
882d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
883999d9058SBarry Smith     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
884b24ad042SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr);
885b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
8867a868f3eSHong Zhang     if (!ijonly){
8877a868f3eSHong Zhang       ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr);
888b24ad042SBarry Smith       ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
889d5b485abSSatish Balay     }
8907a868f3eSHong Zhang   }
891606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
892606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
893d5b485abSSatish Balay 
894d5b485abSSatish Balay   /* Wait on sends1 and sends2 */
895b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr);
896b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr);
897d5b485abSSatish Balay 
8980c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
8990c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
900606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
901606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
902606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
903606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
904d5b485abSSatish Balay 
905d5b485abSSatish Balay   /* Now allocate buffers for a->j, and send them off */
906b24ad042SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr);
907d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
908b24ad042SBarry Smith   ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr);
909d5b485abSSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
910d5b485abSSatish Balay 
911b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr);
912d5b485abSSatish Balay   {
913d5b485abSSatish Balay      for (i=0; i<nrqr; i++) {
914d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
915d5b485abSSatish Balay       sbuf_aj_i = sbuf_aj[i];
916d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0] + 1;
917d5b485abSSatish Balay       ct2       = 0;
918d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
919d5b485abSSatish Balay         kmax = rbuf1[i][2*j];
920d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
921e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
922d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
923d5b485abSSatish Balay           ncols  = nzA + nzB;
924d5b485abSSatish Balay           cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
925d5b485abSSatish Balay 
926d5b485abSSatish Balay           /* load the column indices for this row into cols*/
927d5b485abSSatish Balay           cols  = sbuf_aj_i + ct2;
928d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
929d5b485abSSatish Balay             if ((ctmp = bmap[cworkB[l]]) < cstart)  cols[l] = ctmp;
930d5b485abSSatish Balay             else break;
931d5b485abSSatish Balay           }
932d5b485abSSatish Balay           imark = l;
933d5b485abSSatish Balay           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
934d5b485abSSatish Balay           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
935d5b485abSSatish Balay           ct2 += ncols;
936d5b485abSSatish Balay         }
937d5b485abSSatish Balay       }
938b24ad042SBarry Smith       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
939d5b485abSSatish Balay     }
940d5b485abSSatish Balay   }
941b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr);
942b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr);
943d5b485abSSatish Balay 
944d5b485abSSatish Balay   /* Allocate buffers for a->a, and send them off */
9457a868f3eSHong Zhang   if (!ijonly){
94682502324SSatish Balay     ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr);
947d5b485abSSatish Balay     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
94882502324SSatish Balay     ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr);
949a2feddc5SSatish Balay     for (i=1; i<nrqr; i++)  sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
950d5b485abSSatish Balay 
951b0a32e0cSBarry Smith     ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr);
952d5b485abSSatish Balay     {
953d5b485abSSatish Balay       for (i=0; i<nrqr; i++) {
954d5b485abSSatish Balay         rbuf1_i   = rbuf1[i];
955d5b485abSSatish Balay         sbuf_aa_i = sbuf_aa[i];
956d5b485abSSatish Balay         ct1       = 2*rbuf1_i[0]+1;
957d5b485abSSatish Balay         ct2       = 0;
958d5b485abSSatish Balay         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
959d5b485abSSatish Balay           kmax = rbuf1_i[2*j];
960d5b485abSSatish Balay           for (k=0; k<kmax; k++,ct1++) {
961e8e0db45SSatish Balay             row    = rbuf1_i[ct1] - rstart;
962d5b485abSSatish Balay             nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
963d5b485abSSatish Balay             ncols  = nzA + nzB;
964d5b485abSSatish Balay             cworkA = a_j + a_i[row];     cworkB = b_j + b_i[row];
965a2feddc5SSatish Balay             vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2;
966d5b485abSSatish Balay 
967d5b485abSSatish Balay             /* load the column values for this row into vals*/
9685b83ace0SSatish Balay             vals  = sbuf_aa_i+ct2*bs2;
969d5b485abSSatish Balay             for (l=0; l<nzB; l++) {
9703a40ed3dSBarry Smith               if ((bmap[cworkB[l]]) < cstart) {
9713eda8832SBarry Smith                 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9723a40ed3dSBarry Smith               }
973d5b485abSSatish Balay               else break;
974d5b485abSSatish Balay             }
975d5b485abSSatish Balay             imark = l;
9763a40ed3dSBarry Smith             for (l=0; l<nzA; l++) {
9773eda8832SBarry Smith               ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9783a40ed3dSBarry Smith             }
9793a40ed3dSBarry Smith             for (l=imark; l<nzB; l++) {
9803eda8832SBarry Smith               ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9813a40ed3dSBarry Smith             }
982d5b485abSSatish Balay             ct2 += ncols;
983d5b485abSSatish Balay           }
984d5b485abSSatish Balay         }
9853eda8832SBarry Smith         ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
986d5b485abSSatish Balay       }
987d5b485abSSatish Balay     }
988b0a32e0cSBarry Smith     ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr);
989b0a32e0cSBarry Smith     ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr);
9907a868f3eSHong Zhang   }
991533163c2SBarry Smith   ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr);
992606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
993d5b485abSSatish Balay 
994d5b485abSSatish Balay   /* Form the matrix */
995307b7a18SHong Zhang   /* create col map: global col of C -> local col of submatrices */
996d5b485abSSatish Balay   {
9975d0c19d7SBarry Smith     const PetscInt *icol_i;
998233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
9998fa52d88SHong Zhang     ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&cmap);CHKERRQ(ierr);
1000ff0794f7SSatish Balay     for (i=0; i<ismax; i++) {
1001307b7a18SHong Zhang       if (!allcolumns[i]){
10028fa52d88SHong Zhang         ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr);
1003307b7a18SHong Zhang         jmax   = ncol[i];
1004307b7a18SHong Zhang         icol_i = icol[i];
10058fa52d88SHong Zhang         cmap_i = cmap[i];
1006307b7a18SHong Zhang         for (j=0; j<jmax; j++) {
10078fa52d88SHong Zhang           ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1);CHKERRQ(ierr);
1008307b7a18SHong Zhang         }
1009307b7a18SHong Zhang       } else {
10108fa52d88SHong Zhang         cmap[i] = PETSC_NULL;
1011307b7a18SHong Zhang       }
1012233c2ff6SSatish Balay     }
1013233c2ff6SSatish Balay #else
101423bdbc58SBarry Smith     ierr    = PetscMalloc(ismax*sizeof(PetscInt*),&cmap);CHKERRQ(ierr);
1015d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
10168fa52d88SHong Zhang       if (!allcolumns[i]){
10178fa52d88SHong Zhang         ierr = PetscMalloc(c->Nbs*sizeof(PetscInt),&cmap[i]);CHKERRQ(ierr);
10188fa52d88SHong Zhang         ierr = PetscMemzero(cmap[i],c->Nbs*sizeof(PetscInt));CHKERRQ(ierr);
1019d5b485abSSatish Balay         jmax   = ncol[i];
1020d5b485abSSatish Balay         icol_i = icol[i];
1021d5b485abSSatish Balay         cmap_i = cmap[i];
1022d5b485abSSatish Balay         for (j=0; j<jmax; j++) {
1023d5b485abSSatish Balay           cmap_i[icol_i[j]] = j+1;
1024d5b485abSSatish Balay         }
10258fa52d88SHong Zhang       } else { /* allcolumns[i] */
10268fa52d88SHong Zhang         cmap[i] = PETSC_NULL;
10278fa52d88SHong Zhang       }
1028d5b485abSSatish Balay     }
1029307b7a18SHong Zhang #endif
1030d5b485abSSatish Balay   }
1031d5b485abSSatish Balay 
1032d5b485abSSatish Balay   /* Create lens which is required for MatCreate... */
1033d5b485abSSatish Balay   for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
1034052f0c41SBarry Smith   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr);
1035b24ad042SBarry Smith   lens[0] = (PetscInt*)(lens + ismax);
1036b24ad042SBarry Smith   ierr    = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr);
1037d5b485abSSatish Balay   for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
1038d5b485abSSatish Balay 
1039d5b485abSSatish Balay   /* Update lens from local data */
1040d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1041d5b485abSSatish Balay     jmax   = nrow[i];
10428fa52d88SHong Zhang     if (!allcolumns[i]) cmap_i = cmap[i];
1043d5b485abSSatish Balay     irow_i = irow[i];
1044d5b485abSSatish Balay     lens_i = lens[i];
1045d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1046*4da72fa9SHong Zhang       if (allrows[i]){
1047*4da72fa9SHong Zhang         row = j;
1048*4da72fa9SHong Zhang       } else {
1049d5b485abSSatish Balay         row  = irow_i[j];
1050*4da72fa9SHong Zhang       }
1051233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1052899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1053233c2ff6SSatish Balay #else
1054d5b485abSSatish Balay       proc = rtable[row];
1055233c2ff6SSatish Balay #endif
1056d5b485abSSatish Balay       if (proc == rank) {
105736f4e84dSSatish Balay         /* Get indices from matA and then from matB */
105836f4e84dSSatish Balay         row    = row - rstart;
105936f4e84dSSatish Balay         nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
106036f4e84dSSatish Balay         cworkA =  a_j + a_i[row]; cworkB = b_j + b_i[row];
1061307b7a18SHong Zhang         if (!allcolumns[i]){
1062233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1063233c2ff6SSatish Balay           for (k=0; k<nzA; k++) {
10648fa52d88SHong Zhang             ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1065233c2ff6SSatish Balay             if (tt) { lens_i[j]++; }
1066233c2ff6SSatish Balay           }
1067233c2ff6SSatish Balay           for (k=0; k<nzB; k++) {
10688fa52d88SHong Zhang             ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1069233c2ff6SSatish Balay             if (tt) { lens_i[j]++; }
1070233c2ff6SSatish Balay           }
1071307b7a18SHong Zhang 
1072233c2ff6SSatish Balay #else
1073ca161407SBarry Smith           for (k=0; k<nzA; k++) {
107436f4e84dSSatish Balay             if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; }
1075ca161407SBarry Smith           }
1076ca161407SBarry Smith           for (k=0; k<nzB; k++) {
107736f4e84dSSatish Balay             if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; }
1078d5b485abSSatish Balay           }
1079233c2ff6SSatish Balay #endif
1080307b7a18SHong Zhang         } else {/* allcolumns */
1081307b7a18SHong Zhang           lens_i[j] = nzA + nzB;
1082307b7a18SHong Zhang         }
1083d5b485abSSatish Balay       }
1084a2feddc5SSatish Balay     }
1085ca161407SBarry Smith   }
1086233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1087d5b485abSSatish Balay   /* Create row map*/
10888fa52d88SHong Zhang   ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rmap);CHKERRQ(ierr);
1089ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
10908fa52d88SHong Zhang     ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr);
1091233c2ff6SSatish Balay   }
1092233c2ff6SSatish Balay #else
1093233c2ff6SSatish Balay   /* Create row map*/
1094052f0c41SBarry Smith   ierr    = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr);
1095b24ad042SBarry Smith   rmap[0] = (PetscInt*)(rmap + ismax);
1096b24ad042SBarry Smith   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr);
1097233c2ff6SSatish Balay   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;}
1098233c2ff6SSatish Balay #endif
1099d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1100d5b485abSSatish Balay     irow_i = irow[i];
1101d5b485abSSatish Balay     jmax   = nrow[i];
1102233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
11038fa52d88SHong Zhang     rmap_i = rmap[i];
1104233c2ff6SSatish Balay     for (j=0; j<jmax; j++) {
1105*4da72fa9SHong Zhang       if (allrows[i]){
1106*4da72fa9SHong Zhang         ierr = PetscTableAdd(rmap_i,j+1,j+1);CHKERRQ(ierr);
1107*4da72fa9SHong Zhang       } else {
11088fa52d88SHong Zhang         ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1);CHKERRQ(ierr);
1109233c2ff6SSatish Balay       }
1110*4da72fa9SHong Zhang     }
1111233c2ff6SSatish Balay #else
1112233c2ff6SSatish Balay     rmap_i = rmap[i];
1113d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1114*4da72fa9SHong Zhang       if (allrows[i]){
1115*4da72fa9SHong Zhang         rmap_i[j] = j;
1116*4da72fa9SHong Zhang       } else {
1117d5b485abSSatish Balay         rmap_i[irow_i[j]] = j;
1118d5b485abSSatish Balay       }
1119*4da72fa9SHong Zhang     }
1120233c2ff6SSatish Balay #endif
1121d5b485abSSatish Balay   }
1122d5b485abSSatish Balay 
1123d5b485abSSatish Balay   /* Update lens from offproc data */
1124d5b485abSSatish Balay   {
1125b24ad042SBarry Smith     PetscInt    *rbuf2_i,*rbuf3_i,*sbuf1_i;
1126b24ad042SBarry Smith     PetscMPIInt ii;
1127d5b485abSSatish Balay 
1128d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1129b24ad042SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr);
1130b24ad042SBarry Smith       idex   = pa[ii];
1131999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1132d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1133d5b485abSSatish Balay       ct1     = 2*jmax+1;
1134d5b485abSSatish Balay       ct2     = 0;
1135b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1136b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
1137d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1138d5b485abSSatish Balay         is_no   = sbuf1_i[2*j-1];
1139d5b485abSSatish Balay         max1    = sbuf1_i[2*j];
1140d5b485abSSatish Balay         lens_i  = lens[is_no];
11418fa52d88SHong Zhang         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1142d5b485abSSatish Balay 	rmap_i = rmap[is_no];
1143d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1144233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
11458fa52d88SHong Zhang 	  ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1146233c2ff6SSatish Balay           row--;
1147e32f2f54SBarry Smith           if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1148233c2ff6SSatish Balay #else
1149d5b485abSSatish Balay           row  = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
1150233c2ff6SSatish Balay #endif
1151d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1152d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1153307b7a18SHong Zhang             if (!allcolumns[is_no]){
1154233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
11558fa52d88SHong Zhang               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
1156233c2ff6SSatish Balay               if (tt) {
1157233c2ff6SSatish Balay                 lens_i[row]++;
1158233c2ff6SSatish Balay               }
1159233c2ff6SSatish Balay #else
1160d5b485abSSatish Balay               if (cmap_i[rbuf3_i[ct2]]) {
1161d5b485abSSatish Balay                 lens_i[row]++;
1162d5b485abSSatish Balay               }
1163233c2ff6SSatish Balay #endif
1164307b7a18SHong Zhang             } else { /* allcolumns */
1165307b7a18SHong Zhang               lens_i[row]++;
1166307b7a18SHong Zhang             }
1167d5b485abSSatish Balay           }
1168d5b485abSSatish Balay         }
1169d5b485abSSatish Balay       }
1170d5b485abSSatish Balay     }
1171d5b485abSSatish Balay   }
1172606d414cSSatish Balay   ierr = PetscFree(r_status3);CHKERRQ(ierr);
1173606d414cSSatish Balay   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
11740c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
1175606d414cSSatish Balay   ierr = PetscFree(s_status3);CHKERRQ(ierr);
1176606d414cSSatish Balay   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1177d5b485abSSatish Balay 
1178d5b485abSSatish Balay   /* Create the submatrices */
1179d5b485abSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
11807a868f3eSHong Zhang     if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet");
1181d5b485abSSatish Balay     /*
1182d5b485abSSatish Balay         Assumes new rows are same length as the old rows, hence bug!
1183d5b485abSSatish Balay     */
1184d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1185df36731bSBarry Smith       mat = (Mat_SeqBAIJ *)(submats[i]->data);
1186e7e72b3dSBarry 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");
1187b24ad042SBarry Smith       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
1188e7e72b3dSBarry Smith       if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1189d5b485abSSatish Balay       /* Initial matrix as if empty */
1190b24ad042SBarry Smith       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr);
1191d5f3da31SBarry Smith       submats[i]->factortype = C->factortype;
1192d5b485abSSatish Balay     }
1193ca161407SBarry Smith   } else {
11947a868f3eSHong Zhang     PetscInt bs_tmp;
11957a868f3eSHong Zhang     if (ijonly){
11967a868f3eSHong Zhang       bs_tmp = 1;
11977a868f3eSHong Zhang     } else {
11987a868f3eSHong Zhang       bs_tmp = bs;
11997a868f3eSHong Zhang     }
1200d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1201f69a0ea3SMatthew Knepley       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
12027a868f3eSHong Zhang       ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr);
12037adad957SLisandro Dalcin       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
12047a868f3eSHong Zhang       ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr);
12057a868f3eSHong Zhang       ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */
1206d5b485abSSatish Balay     }
1207d5b485abSSatish Balay   }
1208d5b485abSSatish Balay 
1209d5b485abSSatish Balay   /* Assemble the matrices */
1210d5b485abSSatish Balay   /* First assemble the local rows */
1211d5b485abSSatish Balay   {
1212b24ad042SBarry Smith     PetscInt  ilen_row,*imat_ilen,*imat_j,*imat_i;
1213f3e30f2cSMatthew G Knepley     MatScalar *imat_a = PETSC_NULL;
1214d5b485abSSatish Balay 
1215d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1216df36731bSBarry Smith       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1217d5b485abSSatish Balay       imat_ilen = mat->ilen;
1218d5b485abSSatish Balay       imat_j    = mat->j;
1219d5b485abSSatish Balay       imat_i    = mat->i;
12207a868f3eSHong Zhang       if (!ijonly) imat_a = mat->a;
12218fa52d88SHong Zhang       if (!allcolumns[i]) cmap_i = cmap[i];
1222d5b485abSSatish Balay       rmap_i = rmap[i];
1223d5b485abSSatish Balay       irow_i    = irow[i];
1224d5b485abSSatish Balay       jmax      = nrow[i];
1225d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1226*4da72fa9SHong Zhang         if (allrows[i]){
1227*4da72fa9SHong Zhang           row = j;
1228*4da72fa9SHong Zhang         } else {
1229d5b485abSSatish Balay           row      = irow_i[j];
1230*4da72fa9SHong Zhang         }
1231233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1232899cda47SBarry Smith 	ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1233233c2ff6SSatish Balay #else
1234d5b485abSSatish Balay 	proc = rtable[row];
1235233c2ff6SSatish Balay #endif
1236d5b485abSSatish Balay         if (proc == rank) {
123736f4e84dSSatish Balay           row      = row - rstart;
123836f4e84dSSatish Balay           nzA      = a_i[row+1] - a_i[row];
123936f4e84dSSatish Balay           nzB      = b_i[row+1] - b_i[row];
124036f4e84dSSatish Balay           cworkA   = a_j + a_i[row];
124136f4e84dSSatish Balay           cworkB   = b_j + b_i[row];
12427a868f3eSHong Zhang           if (!ijonly){
124336f4e84dSSatish Balay             vworkA = a_a + a_i[row]*bs2;
124436f4e84dSSatish Balay             vworkB = b_a + b_i[row]*bs2;
12457a868f3eSHong Zhang           }
1246233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
12478fa52d88SHong Zhang 	  ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr);
1248233c2ff6SSatish Balay           row--;
1249e32f2f54SBarry Smith           if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1250233c2ff6SSatish Balay #else
125136f4e84dSSatish Balay           row      = rmap_i[row + rstart];
1252233c2ff6SSatish Balay #endif
1253df36731bSBarry Smith           mat_i    = imat_i[row];
12547a868f3eSHong Zhang           if (!ijonly) mat_a = imat_a + mat_i*bs2;
1255d5b485abSSatish Balay           mat_j    = imat_j + mat_i;
125636f4e84dSSatish Balay           ilen_row = imat_ilen[row];
125736f4e84dSSatish Balay 
125836f4e84dSSatish Balay           /* load the column indices for this row into cols*/
1259307b7a18SHong Zhang           if (!allcolumns[i]){
126036f4e84dSSatish Balay             for (l=0; l<nzB; l++) {
126136f4e84dSSatish Balay               if ((ctmp = bmap[cworkB[l]]) < cstart) {
1262233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
12638fa52d88SHong Zhang                 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr);
1264233c2ff6SSatish Balay                 if (tcol) {
1265233c2ff6SSatish Balay #else
126636f4e84dSSatish Balay                 if ((tcol = cmap_i[ctmp])) {
1267233c2ff6SSatish Balay #endif
1268df36731bSBarry Smith                   *mat_j++ = tcol - 1;
12693eda8832SBarry Smith                   ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1270549d3d68SSatish Balay                   mat_a   += bs2;
1271d5b485abSSatish Balay                   ilen_row++;
1272d5b485abSSatish Balay                 }
1273ca161407SBarry Smith               } else break;
127436f4e84dSSatish Balay             }
127536f4e84dSSatish Balay             imark = l;
127636f4e84dSSatish Balay             for (l=0; l<nzA; l++) {
1277233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
12788fa52d88SHong Zhang               ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1279233c2ff6SSatish Balay               if (tcol) {
1280233c2ff6SSatish Balay #else
128136f4e84dSSatish Balay               if ((tcol = cmap_i[cstart + cworkA[l]])) {
1282233c2ff6SSatish Balay #endif
128336f4e84dSSatish Balay                 *mat_j++ = tcol - 1;
12847a868f3eSHong Zhang                 if (!ijonly){
12853eda8832SBarry Smith                   ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1286549d3d68SSatish Balay                   mat_a += bs2;
12877a868f3eSHong Zhang                 }
128836f4e84dSSatish Balay                 ilen_row++;
128936f4e84dSSatish Balay               }
129036f4e84dSSatish Balay             }
129136f4e84dSSatish Balay             for (l=imark; l<nzB; l++) {
1292233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
12938fa52d88SHong Zhang               ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1294233c2ff6SSatish Balay               if (tcol) {
1295233c2ff6SSatish Balay #else
129636f4e84dSSatish Balay               if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1297233c2ff6SSatish Balay #endif
129836f4e84dSSatish Balay                 *mat_j++ = tcol - 1;
12997a868f3eSHong Zhang                 if (!ijonly){
13003eda8832SBarry Smith                   ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1301549d3d68SSatish Balay                   mat_a   += bs2;
13027a868f3eSHong Zhang                 }
130336f4e84dSSatish Balay                 ilen_row++;
130436f4e84dSSatish Balay               }
130536f4e84dSSatish Balay             }
1306307b7a18SHong Zhang           } else { /* allcolumns */
1307307b7a18SHong Zhang             for (l=0; l<nzB; l++) {
1308307b7a18SHong Zhang               if ((ctmp = bmap[cworkB[l]]) < cstart) {
1309307b7a18SHong Zhang                 *mat_j++ = ctmp;
1310307b7a18SHong Zhang                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1311307b7a18SHong Zhang                 mat_a   += bs2;
1312307b7a18SHong Zhang                 ilen_row++;
1313307b7a18SHong Zhang               } else break;
1314307b7a18SHong Zhang             }
1315307b7a18SHong Zhang             imark = l;
1316307b7a18SHong Zhang             for (l=0; l<nzA; l++) {
1317307b7a18SHong Zhang               *mat_j++ = cstart+cworkA[l];
1318307b7a18SHong Zhang               if (!ijonly){
1319307b7a18SHong Zhang                 ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1320307b7a18SHong Zhang                 mat_a += bs2;
1321307b7a18SHong Zhang               }
1322307b7a18SHong Zhang               ilen_row++;
1323307b7a18SHong Zhang             }
1324307b7a18SHong Zhang             for (l=imark; l<nzB; l++) {
1325307b7a18SHong Zhang               *mat_j++ = bmap[cworkB[l]];
1326307b7a18SHong Zhang               if (!ijonly){
1327307b7a18SHong Zhang                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1328307b7a18SHong Zhang                 mat_a   += bs2;
1329307b7a18SHong Zhang               }
1330307b7a18SHong Zhang               ilen_row++;
1331307b7a18SHong Zhang             }
1332307b7a18SHong Zhang           }
1333d5b485abSSatish Balay           imat_ilen[row] = ilen_row;
1334d5b485abSSatish Balay         }
1335d5b485abSSatish Balay       }
1336d5b485abSSatish Balay     }
1337d5b485abSSatish Balay   }
1338d5b485abSSatish Balay 
1339d5b485abSSatish Balay   /*   Now assemble the off proc rows*/
1340d5b485abSSatish Balay   {
1341b24ad042SBarry Smith     PetscInt    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1342b24ad042SBarry Smith     PetscInt    *imat_j,*imat_i;
1343f3e30f2cSMatthew G Knepley     MatScalar   *imat_a = PETSC_NULL,*rbuf4_i = PETSC_NULL;
1344b24ad042SBarry Smith     PetscMPIInt ii;
1345d5b485abSSatish Balay 
1346d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
13477a868f3eSHong Zhang       if (ijonly){
13487a868f3eSHong Zhang         ii = tmp2;
13497a868f3eSHong Zhang       } else {
1350b24ad042SBarry Smith         ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr);
13517a868f3eSHong Zhang       }
1352b24ad042SBarry Smith       idex   = pa[ii];
1353999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1354d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1355d5b485abSSatish Balay       ct1     = 2*jmax + 1;
1356d5b485abSSatish Balay       ct2     = 0;
1357b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1358b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
13597a868f3eSHong Zhang       if (!ijonly) rbuf4_i = rbuf4[ii];
1360d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1361d5b485abSSatish Balay         is_no     = sbuf1_i[2*j-1];
13628fa52d88SHong Zhang         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1363d5b485abSSatish Balay 	rmap_i    = rmap[is_no];
1364df36731bSBarry Smith         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1365d5b485abSSatish Balay         imat_ilen = mat->ilen;
1366d5b485abSSatish Balay         imat_j    = mat->j;
1367d5b485abSSatish Balay         imat_i    = mat->i;
13687a868f3eSHong Zhang         if (!ijonly) imat_a = mat->a;
1369d5b485abSSatish Balay         max1      = sbuf1_i[2*j];
1370d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1371d5b485abSSatish Balay           row   = sbuf1_i[ct1];
1372233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
13738fa52d88SHong Zhang 	  ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr);
1374233c2ff6SSatish Balay           row--;
1375e32f2f54SBarry Smith           if(row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); }
1376233c2ff6SSatish Balay #else
1377d5b485abSSatish Balay           row   = rmap_i[row];
1378233c2ff6SSatish Balay #endif
1379d5b485abSSatish Balay           ilen  = imat_ilen[row];
1380df36731bSBarry Smith           mat_i = imat_i[row];
13817a868f3eSHong Zhang           if (!ijonly) mat_a = imat_a + mat_i*bs2;
1382d5b485abSSatish Balay           mat_j = imat_j + mat_i;
1383d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1384307b7a18SHong Zhang 
1385307b7a18SHong Zhang           if (!allcolumns[is_no]){
1386d5b485abSSatish Balay             for (l=0; l<max2; l++,ct2++) {
1387233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
13888fa52d88SHong Zhang               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1389233c2ff6SSatish Balay               if (tcol) {
1390233c2ff6SSatish Balay #else
1391d5b485abSSatish Balay               if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1392233c2ff6SSatish Balay #endif
1393df36731bSBarry Smith                 *mat_j++    = tcol - 1;
13947a868f3eSHong Zhang                 if (!ijonly){
13953eda8832SBarry Smith                   ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1396549d3d68SSatish Balay                   mat_a += bs2;
13977a868f3eSHong Zhang                 }
1398d5b485abSSatish Balay                 ilen++;
1399d5b485abSSatish Balay               }
1400d5b485abSSatish Balay             }
1401307b7a18SHong Zhang           } else { /* allcolumns */
1402307b7a18SHong Zhang             for (l=0; l<max2; l++,ct2++) {
1403307b7a18SHong Zhang               *mat_j++    = rbuf3_i[ct2];
1404307b7a18SHong Zhang               if (!ijonly){
1405307b7a18SHong Zhang                 ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1406307b7a18SHong Zhang                 mat_a += bs2;
1407307b7a18SHong Zhang               }
1408307b7a18SHong Zhang               ilen++;
1409307b7a18SHong Zhang             }
1410307b7a18SHong Zhang           }
1411d5b485abSSatish Balay           imat_ilen[row] = ilen;
1412d5b485abSSatish Balay         }
1413d5b485abSSatish Balay       }
1414d5b485abSSatish Balay     }
1415d5b485abSSatish Balay   }
14167a868f3eSHong Zhang   if (!ijonly){
1417606d414cSSatish Balay     ierr = PetscFree(r_status4);CHKERRQ(ierr);
1418606d414cSSatish Balay     ierr = PetscFree(r_waits4);CHKERRQ(ierr);
14190c468ba9SBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
1420606d414cSSatish Balay     ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1421606d414cSSatish Balay     ierr = PetscFree(s_status4);CHKERRQ(ierr);
14227a868f3eSHong Zhang   }
1423d5b485abSSatish Balay 
1424d5b485abSSatish Balay   /* Restore the indices */
1425d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1426*4da72fa9SHong Zhang     if (!allrows[i]){
1427d5b485abSSatish Balay       ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1428*4da72fa9SHong Zhang     }
1429307b7a18SHong Zhang     if (!allcolumns[i]){
1430d5b485abSSatish Balay       ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1431d5b485abSSatish Balay     }
1432307b7a18SHong Zhang   }
1433d5b485abSSatish Balay 
1434d5b485abSSatish Balay   /* Destroy allocated memory */
143523bdbc58SBarry Smith #if defined(PETSC_USE_CTABLE)
143623bdbc58SBarry Smith   ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr);
143723bdbc58SBarry Smith #else
143823bdbc58SBarry Smith   ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr);
143923bdbc58SBarry Smith #endif
144023bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
1441606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
1442d5b485abSSatish Balay 
144323bdbc58SBarry Smith   ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
1444606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1445606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1446d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
1447606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1448d5b485abSSatish Balay   }
1449d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1450606d414cSSatish Balay     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1451d5b485abSSatish Balay   }
145223bdbc58SBarry Smith   ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr);
1453606d414cSSatish Balay   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1454606d414cSSatish Balay   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1455606d414cSSatish Balay   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
14567a868f3eSHong Zhang   if (!ijonly) {
14577a868f3eSHong Zhang     for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);}
14587a868f3eSHong Zhang     ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1459606d414cSSatish Balay     ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1460606d414cSSatish Balay     ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
14617a868f3eSHong Zhang   }
1462d5b485abSSatish Balay 
1463233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1464ff0794f7SSatish Balay   for (i=0; i<ismax; i++) {
14658fa52d88SHong Zhang     ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr);
1466233c2ff6SSatish Balay   }
1467233c2ff6SSatish Balay #endif
14688fa52d88SHong Zhang   ierr = PetscFree(rmap);CHKERRQ(ierr);
14698fa52d88SHong Zhang 
14708fa52d88SHong Zhang   for (i=0; i<ismax; i++){
14718fa52d88SHong Zhang     if (!allcolumns[i]){
14728fa52d88SHong Zhang #if defined (PETSC_USE_CTABLE)
14738fa52d88SHong Zhang       ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr);
14748fa52d88SHong Zhang #else
14758fa52d88SHong Zhang       ierr = PetscFree(cmap[i]);CHKERRQ(ierr);
14768fa52d88SHong Zhang #endif
14778fa52d88SHong Zhang     }
14788fa52d88SHong Zhang   }
14798fa52d88SHong Zhang   ierr = PetscFree(cmap);CHKERRQ(ierr);
1480606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1481d5b485abSSatish Balay 
1482d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
148336f4e84dSSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148436f4e84dSSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1485d5b485abSSatish Balay   }
14867a868f3eSHong Zhang 
14877a868f3eSHong Zhang   c->ijonly = PETSC_FALSE; /* set back to the default */
14883a40ed3dSBarry Smith   PetscFunctionReturn(0);
1489d5b485abSSatish Balay }
1490ca161407SBarry Smith 
1491