xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 30bc6a77c816dcf2db4d7202f4edbdbe9fcef653)
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 
14b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov)
15d5b485abSSatish Balay {
166849ba73SBarry Smith   PetscErrorCode ierr;
17d0f46423SBarry Smith   PetscInt       i,N=C->cmap->N, bs=C->rmap->bs;
1836f4e84dSSatish Balay   IS             *is_new;
1936f4e84dSSatish Balay 
203a40ed3dSBarry Smith   PetscFunctionBegin;
21785e854fSJed Brown   ierr = PetscMalloc1(imax,&is_new);CHKERRQ(ierr);
22df36731bSBarry Smith   /* Convert the indices into block format */
2305d8c843SHong Zhang   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr);
2426fbe8dcSKarl Rupp   if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");
25d5b485abSSatish Balay   for (i=0; i<ov; ++i) {
2636f4e84dSSatish Balay     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr);
27d5b485abSSatish Balay   }
286bf464f9SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);}
2905d8c843SHong Zhang   ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr);
306bf464f9SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);}
31606d414cSSatish Balay   ierr = PetscFree(is_new);CHKERRQ(ierr);
323a40ed3dSBarry Smith   PetscFunctionReturn(0);
33d5b485abSSatish Balay }
34d5b485abSSatish Balay 
35d5b485abSSatish Balay /*
36d5b485abSSatish Balay   Sample message format:
37d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
380e9b0e7eSHong Zhang   to index sets is[1], is[5]
39d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
40d5b485abSSatish Balay   -----------
41d5b485abSSatish Balay   mesg [1] = 1 => is[1]
42d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
43d5b485abSSatish Balay   -----------
44d5b485abSSatish Balay   mesg [5] = 5  => is[5]
45d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
46d5b485abSSatish Balay   -----------
47d5b485abSSatish Balay   mesg [7]
480e9b0e7eSHong Zhang   mesg [n]  data(is[1])
49d5b485abSSatish Balay   -----------
50d5b485abSSatish Balay   mesg[n+1]
51d5b485abSSatish Balay   mesg[m]  data(is[5])
52d5b485abSSatish Balay   -----------
53d5b485abSSatish Balay 
54d5b485abSSatish Balay   Notes:
55d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
56d5b485abSSatish Balay   nrqr - no of requests recieved (which have to be or which have been processed
57d5b485abSSatish Balay */
58db41cccfSHong Zhang PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[])
59d5b485abSSatish Balay {
60df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
615d0c19d7SBarry Smith   const PetscInt **idx,*idx_i;
6224c049a4SHong Zhang   PetscInt       *n,*w3,*w4,**data,len;
636849ba73SBarry Smith   PetscErrorCode ierr;
64b24ad042SBarry Smith   PetscMPIInt    size,rank,tag1,tag2,*w2,*w1,nrqr;
65245d216aSHong Zhang   PetscInt       Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr;
668f9f447aSHong Zhang   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p;
67b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1,*onodes2,*olengths2;
68f1af5d2fSBarry Smith   PetscBT        *table;
69d5b485abSSatish Balay   MPI_Comm       comm;
70d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
71d5b485abSSatish Balay   MPI_Status     *s_status,*recv_status;
728f9f447aSHong Zhang   char           *t_p;
73d5b485abSSatish Balay 
743a40ed3dSBarry Smith   PetscFunctionBegin;
75ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
76d5b485abSSatish Balay   size = c->size;
77d5b485abSSatish Balay   rank = c->rank;
78df36731bSBarry Smith   Mbs  = c->Mbs;
79d5b485abSSatish Balay 
80c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
81c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
82c7dd2924SSatish Balay 
83dcca6d9dSJed Brown   ierr = PetscMalloc2(imax+1,&idx,imax,&n);CHKERRQ(ierr);
8424c049a4SHong Zhang 
85d5b485abSSatish Balay   for (i=0; i<imax; i++) {
86d5b485abSSatish Balay     ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr);
87b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
88d5b485abSSatish Balay   }
89d5b485abSSatish Balay 
90d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
91d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
92884e879aSBarry Smith   ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr);
93d5b485abSSatish Balay   for (i=0; i<imax; i++) {
94b24ad042SBarry Smith     ierr  = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
95d5b485abSSatish Balay     idx_i = idx[i];
96d5b485abSSatish Balay     len   = n[i];
97d5b485abSSatish Balay     for (j=0; j<len; j++) {
98d5b485abSSatish Balay       row = idx_i[j];
99f23aa3ddSBarry Smith       if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries");
100ca31476aSJed Brown       ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr);
101d5b485abSSatish Balay       w4[proc]++;
102d5b485abSSatish Balay     }
103d5b485abSSatish Balay     for (j=0; j<size; j++) {
104d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
105d5b485abSSatish Balay     }
106d5b485abSSatish Balay   }
107d5b485abSSatish Balay 
108d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
109d5b485abSSatish Balay   msz      = 0;              /* total mesg length (for all proc */
1100e9b0e7eSHong Zhang   w1[rank] = 0;              /* no mesg sent to itself */
111d5b485abSSatish Balay   w3[rank] = 0;
112d5b485abSSatish Balay   for (i=0; i<size; i++) {
113d5b485abSSatish Balay     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
114d5b485abSSatish Balay   }
115d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
116854ce69bSBarry Smith   ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr);
117d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
118d5b485abSSatish Balay     if (w1[i]) {pa[j] = i; j++;}
119d5b485abSSatish Balay   }
120d5b485abSSatish Balay 
121d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
122d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
123d5b485abSSatish Balay     j      = pa[i];
124d5b485abSSatish Balay     w1[j] += w2[j] + 2*w3[j];
125d5b485abSSatish Balay     msz   += w1[j];
126d5b485abSSatish Balay   }
127d5b485abSSatish Balay 
128c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
129a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
130c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
131d5b485abSSatish Balay 
132c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
133c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr);
134d5b485abSSatish Balay 
135d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
136dcca6d9dSJed Brown   ierr = PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr);CHKERRQ(ierr);
13723bdbc58SBarry Smith   ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr);
13823bdbc58SBarry Smith   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
139d5b485abSSatish Balay   {
140b24ad042SBarry Smith     PetscInt *iptr = tmp,ict  = 0;
141d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
142d5b485abSSatish Balay       j         = pa[i];
143d5b485abSSatish Balay       iptr     +=  ict;
144d5b485abSSatish Balay       outdat[j] = iptr;
145d5b485abSSatish Balay       ict       = w1[j];
146d5b485abSSatish Balay     }
147d5b485abSSatish Balay   }
148d5b485abSSatish Balay 
149d5b485abSSatish Balay   /* Form the outgoing messages */
150d5b485abSSatish Balay   /*plug in the headers*/
151d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
152d5b485abSSatish Balay     j            = pa[i];
153d5b485abSSatish Balay     outdat[j][0] = 0;
154b24ad042SBarry Smith     ierr         = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
155d5b485abSSatish Balay     ptr[j]       = outdat[j] + 2*w3[j] + 1;
156d5b485abSSatish Balay   }
157d5b485abSSatish Balay 
158d5b485abSSatish Balay   /* Memory for doing local proc's work*/
159d5b485abSSatish Balay   {
1601795a4d1SJed Brown     ierr = PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mbs*imax,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,&t_p);CHKERRQ(ierr);
161d5b485abSSatish Balay 
162bbd702dbSSatish Balay     for (i=0; i<imax; i++) {
163b6410449SSatish Balay       table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
164bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
165d5b485abSSatish Balay     }
166d5b485abSSatish Balay   }
167d5b485abSSatish Balay 
168d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
169d5b485abSSatish Balay   {
170b24ad042SBarry Smith     PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j;
171f1af5d2fSBarry Smith     PetscBT  table_i;
172d5b485abSSatish Balay 
173d5b485abSSatish Balay     for (i=0; i<imax; i++) {
174b24ad042SBarry Smith       ierr    = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
175d5b485abSSatish Balay       n_i     = n[i];
176d5b485abSSatish Balay       table_i = table[i];
177d5b485abSSatish Balay       idx_i   = idx[i];
178d5b485abSSatish Balay       data_i  = data[i];
179d5b485abSSatish Balay       isz_i   = isz[i];
180d5b485abSSatish Balay       for (j=0; j<n_i; j++) {   /* parse the indices of each IS */
181d5b485abSSatish Balay         row  = idx_i[j];
182ca31476aSJed Brown         ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr);
183d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
184d5b485abSSatish Balay           ctr[proc]++;
185d5b485abSSatish Balay           *ptr[proc] = row;
186d5b485abSSatish Balay           ptr[proc]++;
187d6b45a43SBarry Smith         } else { /* Update the local table */
18826fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row;
189d5b485abSSatish Balay         }
190d5b485abSSatish Balay       }
191d5b485abSSatish Balay       /* Update the headers for the current IS */
192d5b485abSSatish Balay       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
193d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
194d5b485abSSatish Balay           outdat_j        = outdat[j];
195d5b485abSSatish Balay           k               = ++outdat_j[0];
196d5b485abSSatish Balay           outdat_j[2*k]   = ctr_j;
197d5b485abSSatish Balay           outdat_j[2*k-1] = i;
198d5b485abSSatish Balay         }
199d5b485abSSatish Balay       }
200d5b485abSSatish Balay       isz[i] = isz_i;
201d5b485abSSatish Balay     }
202d5b485abSSatish Balay   }
203d5b485abSSatish Balay 
204d5b485abSSatish Balay   /*  Now  post the sends */
205854ce69bSBarry Smith   ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr);
206d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
207d5b485abSSatish Balay     j    = pa[i];
208b24ad042SBarry Smith     ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr);
209d5b485abSSatish Balay   }
210d5b485abSSatish Balay 
211d5b485abSSatish Balay   /* No longer need the original indices*/
212d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
213d5b485abSSatish Balay     ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr);
214d5b485abSSatish Balay   }
21524c049a4SHong Zhang   ierr = PetscFree2(idx,n);CHKERRQ(ierr);
216d5b485abSSatish Balay 
217d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
2186bf464f9SBarry Smith     ierr = ISDestroy(&is[i]);CHKERRQ(ierr);
219d5b485abSSatish Balay   }
220d5b485abSSatish Balay 
221d5b485abSSatish Balay   /* Do Local work*/
222df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
223d5b485abSSatish Balay 
224d5b485abSSatish Balay   /* Receive messages*/
225854ce69bSBarry Smith   ierr = PetscMalloc1(nrqr+1,&recv_status);CHKERRQ(ierr);
2260c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);}
227d5b485abSSatish Balay 
228854ce69bSBarry Smith   ierr = PetscMalloc1(nrqs+1,&s_status);CHKERRQ(ierr);
2290c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);}
230d5b485abSSatish Balay 
231d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
23223bdbc58SBarry Smith   ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr);
23323bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
234d5b485abSSatish Balay 
235854ce69bSBarry Smith   ierr = PetscMalloc1(nrqr+1,&xdata);CHKERRQ(ierr);
236854ce69bSBarry Smith   ierr = PetscMalloc1(nrqr+1,&isz1);CHKERRQ(ierr);
237df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
238c05d87d6SBarry Smith   ierr = PetscFree(rbuf[0]);CHKERRQ(ierr);
239606d414cSSatish Balay   ierr = PetscFree(rbuf);CHKERRQ(ierr);
240d5b485abSSatish Balay 
241d5b485abSSatish Balay   /* Send the data back*/
242d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
243d5b485abSSatish Balay   {
244b24ad042SBarry Smith     PetscMPIInt *rw1;
245d5b485abSSatish Balay 
246884e879aSBarry Smith     ierr = PetscCalloc1(size,&rw1);CHKERRQ(ierr);
247c7dd2924SSatish Balay 
248d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
249d5b485abSSatish Balay       proc = recv_status[i].MPI_SOURCE;
250e32f2f54SBarry Smith       if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch");
251d5b485abSSatish Balay       rw1[proc] = isz1[i];
252d5b485abSSatish Balay     }
253d5b485abSSatish Balay 
254c7dd2924SSatish Balay     ierr = PetscFree(onodes1);CHKERRQ(ierr);
255c7dd2924SSatish Balay     ierr = PetscFree(olengths1);CHKERRQ(ierr);
256d5b485abSSatish Balay 
257c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
258c7dd2924SSatish Balay     ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr);
25903512dc5SSatish Balay     ierr = PetscFree(rw1);CHKERRQ(ierr);
260c7dd2924SSatish Balay   }
261c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
262c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr);
263d5b485abSSatish Balay 
264d5b485abSSatish Balay   /*  Now  post the sends */
265854ce69bSBarry Smith   ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr);
266d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
267d5b485abSSatish Balay     j    = recv_status[i].MPI_SOURCE;
268b24ad042SBarry Smith     ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr);
269d5b485abSSatish Balay   }
270d5b485abSSatish Balay 
271d5b485abSSatish Balay   /* receive work done on other processors*/
272d5b485abSSatish Balay   {
273b24ad042SBarry Smith     PetscMPIInt idex;
274b24ad042SBarry Smith     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
275f1af5d2fSBarry Smith     PetscBT     table_i;
276d5b485abSSatish Balay     MPI_Status  *status2;
277d5b485abSSatish Balay 
278854ce69bSBarry Smith     ierr = PetscMalloc1(PetscMax(nrqr,nrqs)+1,&status2);CHKERRQ(ierr);
279d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
280999d9058SBarry Smith       ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr);
281d5b485abSSatish Balay       /* Process the message*/
282999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
283d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
284999d9058SBarry Smith       jmax    = rbuf2[idex][0];
285d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
286d5b485abSSatish Balay         max     = rbuf2_i[2*j];
287d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
288d5b485abSSatish Balay         isz_i   = isz[is_no];
289d5b485abSSatish Balay         data_i  = data[is_no];
290d5b485abSSatish Balay         table_i = table[is_no];
291d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
292d5b485abSSatish Balay           row = rbuf2_i[ct1];
29326fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row;
294d5b485abSSatish Balay         }
295d5b485abSSatish Balay         isz[is_no] = isz_i;
296d5b485abSSatish Balay       }
297d5b485abSSatish Balay     }
2980c468ba9SBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);}
299606d414cSSatish Balay     ierr = PetscFree(status2);CHKERRQ(ierr);
300d5b485abSSatish Balay   }
301d5b485abSSatish Balay 
302d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
30370b3c8c7SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
304d5b485abSSatish Balay   }
305d5b485abSSatish Balay 
306c7dd2924SSatish Balay 
307c7dd2924SSatish Balay   ierr = PetscFree(onodes2);CHKERRQ(ierr);
308c7dd2924SSatish Balay   ierr = PetscFree(olengths2);CHKERRQ(ierr);
309c7dd2924SSatish Balay 
310606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
311c05d87d6SBarry Smith   ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr);
312606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
313606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
314606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
315606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
316606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
3178f9f447aSHong Zhang   ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr);
318606d414cSSatish Balay   ierr = PetscFree(s_status);CHKERRQ(ierr);
319606d414cSSatish Balay   ierr = PetscFree(recv_status);CHKERRQ(ierr);
320606d414cSSatish Balay   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
321606d414cSSatish Balay   ierr = PetscFree(xdata);CHKERRQ(ierr);
322606d414cSSatish Balay   ierr = PetscFree(isz1);CHKERRQ(ierr);
3233a40ed3dSBarry Smith   PetscFunctionReturn(0);
324d5b485abSSatish Balay }
325d5b485abSSatish Balay 
326d5b485abSSatish Balay /*
327df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
328d5b485abSSatish Balay        the work on the local processor.
329d5b485abSSatish Balay 
330d5b485abSSatish Balay      Inputs:
331df36731bSBarry Smith       C      - MAT_MPIBAIJ;
332d5b485abSSatish Balay       imax - total no of index sets processed at a time;
333df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
334d5b485abSSatish Balay 
335d5b485abSSatish Balay      Output:
3360e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
337d5b485abSSatish Balay                to each index set;
338d5b485abSSatish Balay       data   - pointer to the solutions
339d5b485abSSatish Balay */
340b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data)
341d5b485abSSatish Balay {
342df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
343d5b485abSSatish Balay   Mat         A  = c->A,B = c->B;
344df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
345b24ad042SBarry Smith   PetscInt    start,end,val,max,rstart,cstart,*ai,*aj;
346b24ad042SBarry Smith   PetscInt    *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
347f1af5d2fSBarry Smith   PetscBT     table_i;
348d5b485abSSatish Balay 
3493a40ed3dSBarry Smith   PetscFunctionBegin;
350899cda47SBarry Smith   rstart = c->rstartbs;
351899cda47SBarry Smith   cstart = c->cstartbs;
352d5b485abSSatish Balay   ai     = a->i;
353df36731bSBarry Smith   aj     = a->j;
354d5b485abSSatish Balay   bi     = b->i;
355df36731bSBarry Smith   bj     = b->j;
356d5b485abSSatish Balay   garray = c->garray;
357d5b485abSSatish Balay 
358d5b485abSSatish Balay 
359d5b485abSSatish Balay   for (i=0; i<imax; i++) {
360d5b485abSSatish Balay     data_i  = data[i];
361d5b485abSSatish Balay     table_i = table[i];
362d5b485abSSatish Balay     isz_i   = isz[i];
363d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
364d5b485abSSatish Balay       row   = data_i[j] - rstart;
365d5b485abSSatish Balay       start = ai[row];
366d5b485abSSatish Balay       end   = ai[row+1];
367d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
368df36731bSBarry Smith         val = aj[k] + cstart;
36926fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val;
370d5b485abSSatish Balay       }
371d5b485abSSatish Balay       start = bi[row];
372d5b485abSSatish Balay       end   = bi[row+1];
373d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
374df36731bSBarry Smith         val = garray[bj[k]];
37526fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val;
376d5b485abSSatish Balay       }
377d5b485abSSatish Balay     }
378d5b485abSSatish Balay     isz[i] = isz_i;
379d5b485abSSatish Balay   }
3803a40ed3dSBarry Smith   PetscFunctionReturn(0);
381d5b485abSSatish Balay }
382d5b485abSSatish Balay /*
383df36731bSBarry Smith       MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages,
384d5b485abSSatish Balay          and return the output
385d5b485abSSatish Balay 
386d5b485abSSatish Balay          Input:
387d5b485abSSatish Balay            C    - the matrix
388d5b485abSSatish Balay            nrqr - no of messages being processed.
389d5b485abSSatish Balay            rbuf - an array of pointers to the recieved requests
390d5b485abSSatish Balay 
391d5b485abSSatish Balay          Output:
392d5b485abSSatish Balay            xdata - array of messages to be sent back
393d5b485abSSatish Balay            isz1  - size of each message
394d5b485abSSatish Balay 
395a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
396d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
3970e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of
398d5b485abSSatish Balay memory is used.
399d5b485abSSatish Balay 
400d5b485abSSatish Balay */
401b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
402d5b485abSSatish Balay {
403df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
404d5b485abSSatish Balay   Mat            A  = c->A,B = c->B;
405df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
4066849ba73SBarry Smith   PetscErrorCode ierr;
407b24ad042SBarry Smith   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
408b24ad042SBarry Smith   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
409d2a212eaSBarry Smith   PetscInt       val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
410b24ad042SBarry Smith   PetscInt       *rbuf_i,kmax,rbuf_0;
411f1af5d2fSBarry Smith   PetscBT        xtable;
412d5b485abSSatish Balay 
4133a40ed3dSBarry Smith   PetscFunctionBegin;
414df36731bSBarry Smith   Mbs    = c->Mbs;
415899cda47SBarry Smith   rstart = c->rstartbs;
416899cda47SBarry Smith   cstart = c->cstartbs;
417d5b485abSSatish Balay   ai     = a->i;
418df36731bSBarry Smith   aj     = a->j;
419d5b485abSSatish Balay   bi     = b->i;
420df36731bSBarry Smith   bj     = b->j;
421d5b485abSSatish Balay   garray = c->garray;
422d5b485abSSatish Balay 
423d5b485abSSatish Balay 
424d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
425d5b485abSSatish Balay     rbuf_i =  rbuf[i];
426d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
427d5b485abSSatish Balay     ct    += rbuf_0;
42826fbe8dcSKarl Rupp     for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j];
429d5b485abSSatish Balay   }
430d5b485abSSatish Balay 
431701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
432701b8814SBarry Smith   else        max1 = 1;
433d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
434785e854fSJed Brown   ierr         = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr);
435d5b485abSSatish Balay   ++no_malloc;
43653b8de81SBarry Smith   ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr);
437b24ad042SBarry Smith   ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr);
438d5b485abSSatish Balay 
439d5b485abSSatish Balay   ct3 = 0;
440d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
441d5b485abSSatish Balay     rbuf_i =  rbuf[i];
442d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
443d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
444d5b485abSSatish Balay     ct2    =  ct1;
445d5b485abSSatish Balay     ct3   += ct1;
446d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
4476831982aSBarry Smith       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
448d5b485abSSatish Balay       oct2 = ct2;
449d5b485abSSatish Balay       kmax = rbuf_i[2*j];
450d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
451d5b485abSSatish Balay         row = rbuf_i[ct1];
452f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
453d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
454b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
455785e854fSJed Brown             ierr         = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr);
456b24ad042SBarry Smith             ierr         = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
457606d414cSSatish Balay             ierr         = PetscFree(xdata[0]);CHKERRQ(ierr);
458d5b485abSSatish Balay             xdata[0]     = tmp;
459d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
46026fbe8dcSKarl Rupp             for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
461d5b485abSSatish Balay           }
462d5b485abSSatish Balay           xdata[i][ct2++] = row;
463d5b485abSSatish Balay           ct3++;
464d5b485abSSatish Balay         }
465d5b485abSSatish Balay       }
466d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
467d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
468d5b485abSSatish Balay         start = ai[row];
469d5b485abSSatish Balay         end   = ai[row+1];
470d5b485abSSatish Balay         for (l=start; l<end; l++) {
471df36731bSBarry Smith           val = aj[l] + cstart;
472f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
473d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
474b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
475785e854fSJed Brown               ierr         = PetscMalloc1(new_estimate,&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;
48026fbe8dcSKarl Rupp               for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
481d5b485abSSatish Balay             }
482d5b485abSSatish Balay             xdata[i][ct2++] = val;
483d5b485abSSatish Balay             ct3++;
484d5b485abSSatish Balay           }
485d5b485abSSatish Balay         }
486d5b485abSSatish Balay         start = bi[row];
487d5b485abSSatish Balay         end   = bi[row+1];
488d5b485abSSatish Balay         for (l=start; l<end; l++) {
489df36731bSBarry Smith           val = garray[bj[l]];
490f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
491d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
492b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
493785e854fSJed Brown               ierr         = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr);
494b24ad042SBarry Smith               ierr         = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
495606d414cSSatish Balay               ierr         = PetscFree(xdata[0]);CHKERRQ(ierr);
496d5b485abSSatish Balay               xdata[0]     = tmp;
497d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
49826fbe8dcSKarl Rupp               for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
499d5b485abSSatish Balay             }
500d5b485abSSatish Balay             xdata[i][ct2++] = val;
501d5b485abSSatish Balay             ct3++;
502d5b485abSSatish Balay           }
503d5b485abSSatish Balay         }
504d5b485abSSatish Balay       }
505d5b485abSSatish Balay       /* Update the header*/
506d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
507d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
508d5b485abSSatish Balay     }
509d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
510d5b485abSSatish Balay     xdata[i+1]  = xdata[i] + ct2;
511d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
512d5b485abSSatish Balay   }
51394bacf5dSBarry Smith   ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr);
5141e2582c4SBarry Smith   ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr);
5153a40ed3dSBarry Smith   PetscFunctionReturn(0);
516d5b485abSSatish Balay }
517d5b485abSSatish Balay 
518b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
519d5b485abSSatish Balay {
52017df9f7cSHong Zhang   IS             *isrow_block,*iscol_block;
521cf4f527aSSatish Balay   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
5226849ba73SBarry Smith   PetscErrorCode ierr;
523121971b2SHong Zhang   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,N=C->cmap->N,bs=C->rmap->bs;
5245dc98d6aSHong Zhang   Mat_SeqBAIJ    *subc;
5255dc98d6aSHong Zhang   Mat_SubMat     *smat;
526a2feddc5SSatish Balay 
5273a40ed3dSBarry Smith   PetscFunctionBegin;
52829dcf524SDmitry Karpeev   /* The compression and expansion should be avoided. Doesn't point
52929dcf524SDmitry Karpeev      out errors, might change the indices, hence buggey */
53017df9f7cSHong Zhang   ierr = PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block);CHKERRQ(ierr);
53117df9f7cSHong Zhang   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block);CHKERRQ(ierr);
53217df9f7cSHong Zhang   ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block);CHKERRQ(ierr);
533cf4f527aSSatish Balay 
534cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
535b24ad042SBarry Smith   nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt));
536329f5518SBarry Smith   if (!nmax) nmax = 1;
5375dc98d6aSHong Zhang 
5385dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
539cf4f527aSSatish Balay     nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0);
540cf4f527aSSatish Balay 
541653e4784SBarry Smith     /* Make sure every processor loops through the nstages */
542b2566f29SBarry Smith     ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr);
5435dc98d6aSHong Zhang 
5445dc98d6aSHong Zhang     ierr = PetscMalloc1(ismax+1,submat);CHKERRQ(ierr);
5455dc98d6aSHong Zhang   } else {
5465dc98d6aSHong Zhang     subc = (Mat_SeqBAIJ*)((*submat)[0]->data);
5475dc98d6aSHong Zhang     smat   = subc->submatis1;
5485dc98d6aSHong Zhang     if (!smat) {
5495dc98d6aSHong Zhang       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"MatGetSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat");
5505dc98d6aSHong Zhang     }
5515dc98d6aSHong Zhang     nstages = smat->nstages;
5525dc98d6aSHong Zhang   }
5535dc98d6aSHong Zhang 
554cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
555cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
556cf4f527aSSatish Balay     else if (pos == ismax) max_no = 0;
557cf4f527aSSatish Balay     else                   max_no = ismax-pos;
558a42ab0d6SHong Zhang 
559c9ffca76SHong Zhang     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,*submat+pos);CHKERRQ(ierr);
560cf4f527aSSatish Balay     pos += max_no;
561cf4f527aSSatish Balay   }
56236f4e84dSSatish Balay 
5635dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX && ismax) {
5645dc98d6aSHong Zhang     /* save nstages for reuse */
5655dc98d6aSHong Zhang     subc = (Mat_SeqBAIJ*)((*submat)[0]->data);
5665dc98d6aSHong Zhang     smat = subc->submatis1;
5675dc98d6aSHong Zhang     if (!smat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"smat does not exit");
5685dc98d6aSHong Zhang     smat->nstages = nstages;
5695dc98d6aSHong Zhang   }
5705dc98d6aSHong Zhang 
57136f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
57217df9f7cSHong Zhang     ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr);
57317df9f7cSHong Zhang     ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr);
57436f4e84dSSatish Balay   }
57517df9f7cSHong Zhang   ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr);
5763a40ed3dSBarry Smith   PetscFunctionReturn(0);
577a2feddc5SSatish Balay }
578a2feddc5SSatish Balay 
579233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
580e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
581233c2ff6SSatish Balay {
582e005ede5SBarry Smith   PetscInt       nGlobalNd = proc_gnode[size];
5834dc2109aSBarry Smith   PetscMPIInt    fproc;
5844dc2109aSBarry Smith   PetscErrorCode ierr;
585233c2ff6SSatish Balay 
586233c2ff6SSatish Balay   PetscFunctionBegin;
5874dc2109aSBarry Smith   ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr);
58823ce1328SBarry Smith   if (fproc > size) fproc = size;
589e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
590e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
591233c2ff6SSatish Balay     else                         fproc++;
592233c2ff6SSatish Balay   }
593e005ede5SBarry Smith   *rank = fproc;
594233c2ff6SSatish Balay   PetscFunctionReturn(0);
595233c2ff6SSatish Balay }
596233c2ff6SSatish Balay #endif
597233c2ff6SSatish Balay 
598a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
599b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
60017df9f7cSHong Zhang 
60117df9f7cSHong Zhang PetscErrorCode MatDestroy_MPIBAIJ_MatGetSubmatrices(Mat C)
60217df9f7cSHong Zhang {
60317df9f7cSHong Zhang   PetscErrorCode ierr;
60417df9f7cSHong Zhang   Mat_SeqBAIJ    *c = (Mat_SeqBAIJ*)C->data;
60517df9f7cSHong Zhang   Mat_SubMat     *submatj = c->submatis1;
60617df9f7cSHong Zhang   PetscInt       i;
60717df9f7cSHong Zhang 
60817df9f7cSHong Zhang   PetscFunctionBegin;
60917df9f7cSHong Zhang   if (!submatj->id) { /* delete data that are linked only to submats[id=0] */
61017df9f7cSHong Zhang     ierr = PetscFree4(submatj->sbuf1,submatj->ptr,submatj->tmp,submatj->ctr);CHKERRQ(ierr);
61117df9f7cSHong Zhang 
61217df9f7cSHong Zhang     for (i=0; i<submatj->nrqr; ++i) {
61317df9f7cSHong Zhang       ierr = PetscFree(submatj->sbuf2[i]);CHKERRQ(ierr);
61417df9f7cSHong Zhang     }
61517df9f7cSHong Zhang     ierr = PetscFree3(submatj->sbuf2,submatj->req_size,submatj->req_source1);CHKERRQ(ierr);
61617df9f7cSHong Zhang 
61717df9f7cSHong Zhang     if (submatj->rbuf1) {
61817df9f7cSHong Zhang       ierr = PetscFree(submatj->rbuf1[0]);CHKERRQ(ierr);
61917df9f7cSHong Zhang       ierr = PetscFree(submatj->rbuf1);CHKERRQ(ierr);
62017df9f7cSHong Zhang     }
62117df9f7cSHong Zhang 
62217df9f7cSHong Zhang     for (i=0; i<submatj->nrqs; ++i) {
62317df9f7cSHong Zhang       ierr = PetscFree(submatj->rbuf3[i]);CHKERRQ(ierr);
62417df9f7cSHong Zhang     }
62517df9f7cSHong Zhang     ierr = PetscFree3(submatj->req_source2,submatj->rbuf2,submatj->rbuf3);CHKERRQ(ierr);
62617df9f7cSHong Zhang     ierr = PetscFree(submatj->pa);CHKERRQ(ierr);
62717df9f7cSHong Zhang   }
62817df9f7cSHong Zhang 
62917df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
63017df9f7cSHong Zhang   ierr = PetscTableDestroy((PetscTable*)&submatj->rmap);CHKERRQ(ierr);
63117df9f7cSHong Zhang   if (submatj->cmap_loc) {ierr = PetscFree(submatj->cmap_loc);CHKERRQ(ierr);}
63217df9f7cSHong Zhang   ierr = PetscFree(submatj->rmap_loc);CHKERRQ(ierr);
63317df9f7cSHong Zhang #else
63417df9f7cSHong Zhang   ierr = PetscFree(submatj->rmap);CHKERRQ(ierr);
63517df9f7cSHong Zhang #endif
63617df9f7cSHong Zhang 
63717df9f7cSHong Zhang   if (!submatj->allcolumns) {
63817df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
63917df9f7cSHong Zhang     ierr = PetscTableDestroy((PetscTable*)&submatj->cmap);CHKERRQ(ierr);
64017df9f7cSHong Zhang #else
64117df9f7cSHong Zhang     ierr = PetscFree(submatj->cmap);CHKERRQ(ierr);
64217df9f7cSHong Zhang #endif
64317df9f7cSHong Zhang   }
64417df9f7cSHong Zhang   ierr = submatj->destroy(C);CHKERRQ(ierr);
64517df9f7cSHong Zhang   ierr = PetscFree(submatj->row2proc);CHKERRQ(ierr);
64617df9f7cSHong Zhang 
64717df9f7cSHong Zhang   ierr = PetscFree(submatj);CHKERRQ(ierr);
64817df9f7cSHong Zhang   PetscFunctionReturn(0);
64917df9f7cSHong Zhang }
65017df9f7cSHong Zhang 
651c9ffca76SHong Zhang PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
65217df9f7cSHong Zhang {
65317df9f7cSHong Zhang   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
65417df9f7cSHong Zhang   Mat            A  = c->A;
65517df9f7cSHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc;
65617df9f7cSHong Zhang   const PetscInt **icol,**irow;
65717df9f7cSHong Zhang   PetscInt       *nrow,*ncol,start;
65817df9f7cSHong Zhang   PetscErrorCode ierr;
65917df9f7cSHong Zhang   PetscMPIInt    rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr;
66017df9f7cSHong Zhang   PetscInt       **sbuf1,**sbuf2,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1;
66117df9f7cSHong Zhang   PetscInt       nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol;
66217df9f7cSHong Zhang   PetscInt       **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2;
66317df9f7cSHong Zhang   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
66417df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
66517df9f7cSHong Zhang   PetscTable     *cmap,cmap_i=NULL,*rmap,rmap_i;
66617df9f7cSHong Zhang #else
66717df9f7cSHong Zhang   PetscInt       **cmap,*cmap_i=NULL,**rmap,*rmap_i;
66817df9f7cSHong Zhang #endif
66917df9f7cSHong Zhang   const PetscInt *irow_i;
67017df9f7cSHong Zhang   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
67117df9f7cSHong Zhang   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
67217df9f7cSHong Zhang   MPI_Request    *r_waits4,*s_waits3,*s_waits4;
67317df9f7cSHong Zhang   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
67417df9f7cSHong Zhang   MPI_Status     *r_status3,*r_status4,*s_status4;
67517df9f7cSHong Zhang   MPI_Comm       comm;
676*30bc6a77SHong Zhang   PetscScalar    **rbuf4,*rbuf4_i=NULL,**sbuf_aa,*vals,*mat_a=NULL,*imat_a,*sbuf_aa_i;
67717df9f7cSHong Zhang   PetscMPIInt    *onodes1,*olengths1,end;
67880d31651SHong Zhang   PetscInt       **row2proc,*row2proc_i,*imat_ilen,*imat_j,*imat_i;
67917df9f7cSHong Zhang   Mat_SubMat     **smats,*smat_i;
68017df9f7cSHong Zhang   PetscBool      *issorted,colflag,iscsorted=PETSC_TRUE;
68117df9f7cSHong Zhang   PetscInt       *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen;
682a42ab0d6SHong Zhang   PetscInt       bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs;
683a42ab0d6SHong Zhang   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
684a42ab0d6SHong Zhang   PetscInt       nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB;
685afb3cbd8SHong Zhang   PetscScalar    *vworkA=NULL,*vworkB=NULL,*a_a = a->a,*b_a = b->a;
686a42ab0d6SHong Zhang   PetscInt       cstart = c->cstartbs,*bmap = c->garray;
68780d31651SHong Zhang   PetscBool      *allrows,*allcolumns;
688a42ab0d6SHong Zhang 
68917df9f7cSHong Zhang   PetscFunctionBegin;
69017df9f7cSHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
69117df9f7cSHong Zhang   size = c->size;
69217df9f7cSHong Zhang   rank = c->rank;
69317df9f7cSHong Zhang 
69480d31651SHong Zhang   ierr = PetscMalloc6(ismax,&smats,ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax,&allcolumns,ismax,&allrows);CHKERRQ(ierr);
69517df9f7cSHong Zhang   ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr);
69617df9f7cSHong Zhang 
69717df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
69817df9f7cSHong Zhang     ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr);
69980d31651SHong Zhang     if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */
70017df9f7cSHong Zhang     ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr);
70117df9f7cSHong Zhang 
702ec3c739cSHong Zhang     /* Check for special case: allcolumns */
70317df9f7cSHong Zhang     ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr);
70417df9f7cSHong Zhang     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
7055dd0c08cSHong Zhang 
7061abc2f7bSHong Zhang     if (colflag && ncol[i] == c->Nbs) {
7071abc2f7bSHong Zhang       allcolumns[i] = PETSC_TRUE;
7081abc2f7bSHong Zhang       icol[i]       = NULL;
7091abc2f7bSHong Zhang     } else {
71017df9f7cSHong Zhang       allcolumns[i] = PETSC_FALSE;
71117df9f7cSHong Zhang       ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
7121abc2f7bSHong Zhang     }
713ec3c739cSHong Zhang 
714ec3c739cSHong Zhang     /* Check for special case: allrows */
715ec3c739cSHong Zhang     ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr);
716ec3c739cSHong Zhang     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
717ec3c739cSHong Zhang     if (colflag && nrow[i] == c->Mbs) {
718ec3c739cSHong Zhang       allrows[i] = PETSC_TRUE;
719ec3c739cSHong Zhang       irow[i]    = NULL;
720ec3c739cSHong Zhang     } else {
721ec3c739cSHong Zhang       allrows[i] = PETSC_FALSE;
722ec3c739cSHong Zhang       ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
723ec3c739cSHong Zhang     }
72417df9f7cSHong Zhang   }
72517df9f7cSHong Zhang 
72617df9f7cSHong Zhang   if (scall == MAT_REUSE_MATRIX) {
72717df9f7cSHong Zhang     /* Assumes new rows are same length as the old rows */
72817df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
72917df9f7cSHong Zhang       subc = (Mat_SeqBAIJ*)(submats[i]->data);
7309e686edcSHong Zhang       if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
73117df9f7cSHong Zhang 
73217df9f7cSHong Zhang       /* Initial matrix as if empty */
7339e686edcSHong Zhang       ierr = PetscMemzero(subc->ilen,subc->mbs*sizeof(PetscInt));CHKERRQ(ierr);
73417df9f7cSHong Zhang 
73517df9f7cSHong Zhang       /* Initial matrix as if empty */
73617df9f7cSHong Zhang       submats[i]->factortype = C->factortype;
73717df9f7cSHong Zhang 
73817df9f7cSHong Zhang       smat_i   = subc->submatis1;
73917df9f7cSHong Zhang       smats[i] = smat_i;
74017df9f7cSHong Zhang 
74117df9f7cSHong Zhang       nrqs        = smat_i->nrqs;
74217df9f7cSHong Zhang       nrqr        = smat_i->nrqr;
74317df9f7cSHong Zhang       rbuf1       = smat_i->rbuf1;
74417df9f7cSHong Zhang       rbuf2       = smat_i->rbuf2;
74517df9f7cSHong Zhang       rbuf3       = smat_i->rbuf3;
74617df9f7cSHong Zhang       req_source2 = smat_i->req_source2;
74717df9f7cSHong Zhang 
74817df9f7cSHong Zhang       sbuf1     = smat_i->sbuf1;
74917df9f7cSHong Zhang       sbuf2     = smat_i->sbuf2;
75017df9f7cSHong Zhang       ptr       = smat_i->ptr;
75117df9f7cSHong Zhang       tmp       = smat_i->tmp;
75217df9f7cSHong Zhang       ctr       = smat_i->ctr;
75317df9f7cSHong Zhang 
75417df9f7cSHong Zhang       pa          = smat_i->pa;
75517df9f7cSHong Zhang       req_size    = smat_i->req_size;
75617df9f7cSHong Zhang       req_source1 = smat_i->req_source1;
75717df9f7cSHong Zhang 
75817df9f7cSHong Zhang       allcolumns[i] = smat_i->allcolumns;
759ec3c739cSHong Zhang       allrows[i]    = smat_i->allrows;
76017df9f7cSHong Zhang       row2proc[i]   = smat_i->row2proc;
76117df9f7cSHong Zhang       rmap[i]       = smat_i->rmap;
76217df9f7cSHong Zhang       cmap[i]       = smat_i->cmap;
76317df9f7cSHong Zhang     }
76417df9f7cSHong Zhang   } else { /* scall == MAT_INITIAL_MATRIX */
76517df9f7cSHong Zhang     /* Get some new tags to keep the communication clean */
76617df9f7cSHong Zhang     ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
76717df9f7cSHong Zhang     ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
76817df9f7cSHong Zhang 
76917df9f7cSHong Zhang     /* evaluate communication - mesg to who, length of mesg, and buffer space
77017df9f7cSHong Zhang      required. Based on this, buffers are allocated, and data copied into them*/
77117df9f7cSHong Zhang     ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr);   /* mesg size, initialize work vectors */
77217df9f7cSHong Zhang 
77317df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
77417df9f7cSHong Zhang       jmax   = nrow[i];
77517df9f7cSHong Zhang       irow_i = irow[i];
77617df9f7cSHong Zhang 
77717df9f7cSHong Zhang       ierr   = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr);
77817df9f7cSHong Zhang       row2proc[i] = row2proc_i;
77917df9f7cSHong Zhang 
78017df9f7cSHong Zhang       if (issorted[i]) proc = 0;
78117df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
78217df9f7cSHong Zhang         if (!issorted[i]) proc = 0;
783ec3c739cSHong Zhang         if (allrows[i]) row = j;
784ec3c739cSHong Zhang         else row = irow_i[j];
785ec3c739cSHong Zhang 
786a42ab0d6SHong Zhang         while (row >= c->rangebs[proc+1]) proc++;
78717df9f7cSHong Zhang         w4[proc]++;
78817df9f7cSHong Zhang         row2proc_i[j] = proc; /* map row index to proc */
78917df9f7cSHong Zhang       }
79017df9f7cSHong Zhang       for (j=0; j<size; j++) {
79117df9f7cSHong Zhang         if (w4[j]) { w1[j] += w4[j];  w3[j]++; w4[j] = 0;}
79217df9f7cSHong Zhang       }
79317df9f7cSHong Zhang     }
79417df9f7cSHong Zhang 
79517df9f7cSHong Zhang     nrqs     = 0;              /* no of outgoing messages */
79617df9f7cSHong Zhang     msz      = 0;              /* total mesg length (for all procs) */
79717df9f7cSHong Zhang     w1[rank] = 0;              /* no mesg sent to self */
79817df9f7cSHong Zhang     w3[rank] = 0;
79917df9f7cSHong Zhang     for (i=0; i<size; i++) {
80017df9f7cSHong Zhang       if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
80117df9f7cSHong Zhang     }
80217df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/
80317df9f7cSHong Zhang     for (i=0,j=0; i<size; i++) {
80417df9f7cSHong Zhang       if (w1[i]) { pa[j] = i; j++; }
80517df9f7cSHong Zhang     }
80617df9f7cSHong Zhang 
80717df9f7cSHong Zhang     /* Each message would have a header = 1 + 2*(no of IS) + data */
80817df9f7cSHong Zhang     for (i=0; i<nrqs; i++) {
80917df9f7cSHong Zhang       j      = pa[i];
81017df9f7cSHong Zhang       w1[j] += w2[j] + 2* w3[j];
81117df9f7cSHong Zhang       msz   += w1[j];
81217df9f7cSHong Zhang     }
81317df9f7cSHong Zhang     ierr = PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);CHKERRQ(ierr);
81417df9f7cSHong Zhang 
81517df9f7cSHong Zhang     /* Determine the number of messages to expect, their lengths, from from-ids */
81617df9f7cSHong Zhang     ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
81717df9f7cSHong Zhang     ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
81817df9f7cSHong Zhang 
81917df9f7cSHong Zhang     /* Now post the Irecvs corresponding to these messages */
82017df9f7cSHong Zhang     tag0 = ((PetscObject)C)->tag;
82117df9f7cSHong Zhang     ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
82217df9f7cSHong Zhang 
82317df9f7cSHong Zhang     ierr = PetscFree(onodes1);CHKERRQ(ierr);
82417df9f7cSHong Zhang     ierr = PetscFree(olengths1);CHKERRQ(ierr);
82517df9f7cSHong Zhang 
82617df9f7cSHong Zhang     /* Allocate Memory for outgoing messages */
82717df9f7cSHong Zhang     ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr);
82817df9f7cSHong Zhang     ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
82917df9f7cSHong Zhang     ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
83017df9f7cSHong Zhang 
83117df9f7cSHong Zhang     {
83217df9f7cSHong Zhang       PetscInt *iptr = tmp;
83317df9f7cSHong Zhang       k    = 0;
83417df9f7cSHong Zhang       for (i=0; i<nrqs; i++) {
83517df9f7cSHong Zhang         j        = pa[i];
83617df9f7cSHong Zhang         iptr    += k;
83717df9f7cSHong Zhang         sbuf1[j] = iptr;
83817df9f7cSHong Zhang         k        = w1[j];
83917df9f7cSHong Zhang       }
84017df9f7cSHong Zhang     }
84117df9f7cSHong Zhang 
84217df9f7cSHong Zhang     /* Form the outgoing messages. Initialize the header space */
84317df9f7cSHong Zhang     for (i=0; i<nrqs; i++) {
84417df9f7cSHong Zhang       j           = pa[i];
84517df9f7cSHong Zhang       sbuf1[j][0] = 0;
84617df9f7cSHong Zhang       ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
84717df9f7cSHong Zhang       ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
84817df9f7cSHong Zhang     }
84917df9f7cSHong Zhang 
85017df9f7cSHong Zhang     /* Parse the isrow and copy data into outbuf */
85117df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
85217df9f7cSHong Zhang       row2proc_i = row2proc[i];
85317df9f7cSHong Zhang       ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
85417df9f7cSHong Zhang       irow_i = irow[i];
85517df9f7cSHong Zhang       jmax   = nrow[i];
85617df9f7cSHong Zhang       for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
85717df9f7cSHong Zhang         proc = row2proc_i[j];
858ec3c739cSHong Zhang         if (allrows[i]) row = j;
859ec3c739cSHong Zhang         else row = irow_i[j];
860ec3c739cSHong Zhang 
86117df9f7cSHong Zhang         if (proc != rank) { /* copy to the outgoing buf*/
86217df9f7cSHong Zhang           ctr[proc]++;
863ec3c739cSHong Zhang           *ptr[proc] = row;
86417df9f7cSHong Zhang           ptr[proc]++;
86517df9f7cSHong Zhang         }
86617df9f7cSHong Zhang       }
86717df9f7cSHong Zhang       /* Update the headers for the current IS */
86817df9f7cSHong Zhang       for (j=0; j<size; j++) { /* Can Optimise this loop too */
86917df9f7cSHong Zhang         if ((ctr_j = ctr[j])) {
87017df9f7cSHong Zhang           sbuf1_j        = sbuf1[j];
87117df9f7cSHong Zhang           k              = ++sbuf1_j[0];
87217df9f7cSHong Zhang           sbuf1_j[2*k]   = ctr_j;
87317df9f7cSHong Zhang           sbuf1_j[2*k-1] = i;
87417df9f7cSHong Zhang         }
87517df9f7cSHong Zhang       }
87617df9f7cSHong Zhang     }
87717df9f7cSHong Zhang 
87817df9f7cSHong Zhang     /*  Now  post the sends */
87917df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr);
88017df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
88117df9f7cSHong Zhang       j    = pa[i];
88217df9f7cSHong Zhang       ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
88317df9f7cSHong Zhang     }
88417df9f7cSHong Zhang 
88517df9f7cSHong Zhang     /* Post Receives to capture the buffer size */
88617df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr);
88717df9f7cSHong Zhang     ierr = PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);CHKERRQ(ierr);
88817df9f7cSHong Zhang     rbuf2[0] = tmp + msz;
88917df9f7cSHong Zhang     for (i=1; i<nrqs; ++i) {
89017df9f7cSHong Zhang       rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
89117df9f7cSHong Zhang     }
89217df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
89317df9f7cSHong Zhang       j    = pa[i];
89417df9f7cSHong Zhang       ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRQ(ierr);
89517df9f7cSHong Zhang     }
89617df9f7cSHong Zhang 
89717df9f7cSHong Zhang     /* Send to other procs the buf size they should allocate */
89817df9f7cSHong Zhang     /* Receive messages*/
89917df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr);
90017df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr);
90117df9f7cSHong Zhang     ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr);
90217df9f7cSHong Zhang     {
903a42ab0d6SHong Zhang       PetscInt   *sAi = a->i,*sBi = b->i,id;
90417df9f7cSHong Zhang       PetscInt   *sbuf2_i;
90517df9f7cSHong Zhang 
90617df9f7cSHong Zhang       ierr = MPI_Waitall(nrqr,r_waits1,r_status1);CHKERRQ(ierr);
90717df9f7cSHong Zhang       for (i=0; i<nrqr; ++i) {
90817df9f7cSHong Zhang         req_size[i] = 0;
90917df9f7cSHong Zhang         rbuf1_i        = rbuf1[i];
91017df9f7cSHong Zhang         start          = 2*rbuf1_i[0] + 1;
91117df9f7cSHong Zhang         ierr           = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
91217df9f7cSHong Zhang         ierr           = PetscMalloc1(end+1,&sbuf2[i]);CHKERRQ(ierr);
91317df9f7cSHong Zhang         sbuf2_i        = sbuf2[i];
91417df9f7cSHong Zhang         for (j=start; j<end; j++) {
91517df9f7cSHong Zhang           id              = rbuf1_i[j] - rstart;
91617df9f7cSHong Zhang           ncols           = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
91717df9f7cSHong Zhang           sbuf2_i[j]      = ncols;
91817df9f7cSHong Zhang           req_size[i] += ncols;
91917df9f7cSHong Zhang         }
92017df9f7cSHong Zhang         req_source1[i] = r_status1[i].MPI_SOURCE;
92117df9f7cSHong Zhang         /* form the header */
92217df9f7cSHong Zhang         sbuf2_i[0] = req_size[i];
92317df9f7cSHong Zhang         for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
92417df9f7cSHong Zhang 
92517df9f7cSHong Zhang         ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRQ(ierr);
92617df9f7cSHong Zhang       }
92717df9f7cSHong Zhang     }
92817df9f7cSHong Zhang     ierr = PetscFree(r_status1);CHKERRQ(ierr);
92917df9f7cSHong Zhang     ierr = PetscFree(r_waits1);CHKERRQ(ierr);
93017df9f7cSHong Zhang     ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
93117df9f7cSHong Zhang 
93217df9f7cSHong Zhang     /* Receive messages*/
93317df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr);
93417df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr);
93517df9f7cSHong Zhang 
93617df9f7cSHong Zhang     ierr = MPI_Waitall(nrqs,r_waits2,r_status2);CHKERRQ(ierr);
93717df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
93817df9f7cSHong Zhang       ierr = PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);CHKERRQ(ierr);
93917df9f7cSHong Zhang       req_source2[i] = r_status2[i].MPI_SOURCE;
94017df9f7cSHong Zhang       ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRQ(ierr);
94117df9f7cSHong Zhang     }
94217df9f7cSHong Zhang     ierr = PetscFree(r_status2);CHKERRQ(ierr);
94317df9f7cSHong Zhang     ierr = PetscFree(r_waits2);CHKERRQ(ierr);
94417df9f7cSHong Zhang 
94517df9f7cSHong Zhang     /* Wait on sends1 and sends2 */
94617df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr);
94717df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr);
94817df9f7cSHong Zhang 
94917df9f7cSHong Zhang     if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
95017df9f7cSHong Zhang     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
95117df9f7cSHong Zhang     ierr = PetscFree(s_status1);CHKERRQ(ierr);
95217df9f7cSHong Zhang     ierr = PetscFree(s_status2);CHKERRQ(ierr);
95317df9f7cSHong Zhang     ierr = PetscFree(s_waits1);CHKERRQ(ierr);
95417df9f7cSHong Zhang     ierr = PetscFree(s_waits2);CHKERRQ(ierr);
95517df9f7cSHong Zhang 
95617df9f7cSHong Zhang     /* Now allocate sending buffers for a->j, and send them off */
95717df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr);
95817df9f7cSHong Zhang     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
95917df9f7cSHong Zhang     ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr);
96017df9f7cSHong Zhang     for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
96117df9f7cSHong Zhang 
96217df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr);
96317df9f7cSHong Zhang     {
96417df9f7cSHong Zhang 
96517df9f7cSHong Zhang       for (i=0; i<nrqr; i++) {
96617df9f7cSHong Zhang         rbuf1_i   = rbuf1[i];
96717df9f7cSHong Zhang         sbuf_aj_i = sbuf_aj[i];
96817df9f7cSHong Zhang         ct1       = 2*rbuf1_i[0] + 1;
96917df9f7cSHong Zhang         ct2       = 0;
97017df9f7cSHong Zhang         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
97117df9f7cSHong Zhang           kmax = rbuf1[i][2*j];
97217df9f7cSHong Zhang           for (k=0; k<kmax; k++,ct1++) {
97317df9f7cSHong Zhang             row    = rbuf1_i[ct1] - rstart;
97417df9f7cSHong Zhang             nzA    = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row];
97517df9f7cSHong Zhang             ncols  = nzA + nzB;
97617df9f7cSHong Zhang             cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
97717df9f7cSHong Zhang 
97817df9f7cSHong Zhang             /* load the column indices for this row into cols */
97917df9f7cSHong Zhang             cols = sbuf_aj_i + ct2;
98017df9f7cSHong Zhang             for (l=0; l<nzB; l++) {
981a42ab0d6SHong Zhang               if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
982a42ab0d6SHong Zhang               else break;
98317df9f7cSHong Zhang             }
984a42ab0d6SHong Zhang             imark = l;
985a42ab0d6SHong Zhang             for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];}
986a42ab0d6SHong Zhang             for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
98717df9f7cSHong Zhang             ct2 += ncols;
98817df9f7cSHong Zhang           }
98917df9f7cSHong Zhang         }
99017df9f7cSHong Zhang         ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRQ(ierr);
99117df9f7cSHong Zhang       }
99217df9f7cSHong Zhang     }
99317df9f7cSHong Zhang     ierr = PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);CHKERRQ(ierr);
99417df9f7cSHong Zhang 
99517df9f7cSHong Zhang     /* create col map: global col of C -> local col of submatrices */
99617df9f7cSHong Zhang     const PetscInt *icol_i;
99717df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
99817df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
99917df9f7cSHong Zhang       if (!allcolumns[i]) {
1000a42ab0d6SHong Zhang         ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr);
100117df9f7cSHong Zhang 
100217df9f7cSHong Zhang         jmax   = ncol[i];
100317df9f7cSHong Zhang         icol_i = icol[i];
100417df9f7cSHong Zhang         cmap_i = cmap[i];
100517df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
100617df9f7cSHong Zhang           ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
100717df9f7cSHong Zhang         }
100817df9f7cSHong Zhang       } else cmap[i] = NULL;
100917df9f7cSHong Zhang     }
101017df9f7cSHong Zhang #else
101117df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
101217df9f7cSHong Zhang       if (!allcolumns[i]) {
1013a42ab0d6SHong Zhang         ierr   = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr);
101417df9f7cSHong Zhang         jmax   = ncol[i];
101517df9f7cSHong Zhang         icol_i = icol[i];
101617df9f7cSHong Zhang         cmap_i = cmap[i];
1017a42ab0d6SHong Zhang         for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1;
101817df9f7cSHong Zhang       } else cmap[i] = NULL;
101917df9f7cSHong Zhang     }
102017df9f7cSHong Zhang #endif
102117df9f7cSHong Zhang 
102217df9f7cSHong Zhang     /* Create lens which is required for MatCreate... */
102317df9f7cSHong Zhang     for (i=0,j=0; i<ismax; i++) j += nrow[i];
102417df9f7cSHong Zhang     ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr);
102517df9f7cSHong Zhang 
102617df9f7cSHong Zhang     if (ismax) {
102717df9f7cSHong Zhang       ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr);
102817df9f7cSHong Zhang     }
102917df9f7cSHong Zhang     for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];
103017df9f7cSHong Zhang 
103117df9f7cSHong Zhang     /* Update lens from local data */
103217df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
103317df9f7cSHong Zhang       row2proc_i = row2proc[i];
103417df9f7cSHong Zhang       jmax = nrow[i];
103517df9f7cSHong Zhang       if (!allcolumns[i]) cmap_i = cmap[i];
103617df9f7cSHong Zhang       irow_i = irow[i];
103717df9f7cSHong Zhang       lens_i = lens[i];
103817df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
1039ec3c739cSHong Zhang         if (allrows[i]) row = j;
1040ec3c739cSHong Zhang         else row = irow_i[j]; /* global blocked row of C */
1041ec3c739cSHong Zhang 
104217df9f7cSHong Zhang         proc = row2proc_i[j];
104317df9f7cSHong Zhang         if (proc == rank) {
1044a42ab0d6SHong Zhang           /* Get indices from matA and then from matB */
104517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1046a42ab0d6SHong Zhang           PetscInt   tt;
1047a42ab0d6SHong Zhang #endif
1048a42ab0d6SHong Zhang           row    = row - rstart;
1049a42ab0d6SHong Zhang           nzA    = a_i[row+1] - a_i[row];
1050a42ab0d6SHong Zhang           nzB    = b_i[row+1] - b_i[row];
1051a42ab0d6SHong Zhang           cworkA =  a_j + a_i[row];
1052a42ab0d6SHong Zhang           cworkB = b_j + b_i[row];
1053a42ab0d6SHong Zhang 
1054a42ab0d6SHong Zhang           if (!allcolumns[i]) {
1055a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1056a42ab0d6SHong Zhang             for (k=0; k<nzA; k++) {
1057a42ab0d6SHong Zhang               ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1058a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1059a42ab0d6SHong Zhang             }
1060a42ab0d6SHong Zhang             for (k=0; k<nzB; k++) {
1061a42ab0d6SHong Zhang               ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1062a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1063a42ab0d6SHong Zhang             }
1064a42ab0d6SHong Zhang 
1065a42ab0d6SHong Zhang #else
1066a42ab0d6SHong Zhang             for (k=0; k<nzA; k++) {
1067a42ab0d6SHong Zhang               if (cmap_i[cstart + cworkA[k]]) lens_i[j]++;
1068a42ab0d6SHong Zhang             }
1069a42ab0d6SHong Zhang             for (k=0; k<nzB; k++) {
1070a42ab0d6SHong Zhang               if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++;
1071a42ab0d6SHong Zhang             }
1072a42ab0d6SHong Zhang #endif
1073a42ab0d6SHong Zhang           } else { /* allcolumns */
1074a42ab0d6SHong Zhang             lens_i[j] = nzA + nzB;
1075a42ab0d6SHong Zhang           }
107617df9f7cSHong Zhang         }
107717df9f7cSHong Zhang       }
107817df9f7cSHong Zhang     }
107917df9f7cSHong Zhang 
108017df9f7cSHong Zhang     /* Create row map: global row of C -> local row of submatrices */
108117df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
1082059f6b74SHong Zhang       if (!allrows[i]) {
1083059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE)
1084a42ab0d6SHong Zhang         ierr   = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr);
108517df9f7cSHong Zhang         irow_i = irow[i];
108617df9f7cSHong Zhang         jmax   = nrow[i];
108717df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
1088ec3c739cSHong Zhang           if (allrows[i]) {
1089ec3c739cSHong Zhang             ierr = PetscTableAdd(rmap[i],j+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1090ec3c739cSHong Zhang           } else {
109117df9f7cSHong Zhang             ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
109217df9f7cSHong Zhang           }
109317df9f7cSHong Zhang         }
109417df9f7cSHong Zhang #else
1095a42ab0d6SHong Zhang         ierr   = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr);
109617df9f7cSHong Zhang         rmap_i = rmap[i];
109717df9f7cSHong Zhang         irow_i = irow[i];
109817df9f7cSHong Zhang         jmax   = nrow[i];
109917df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
1100ec3c739cSHong Zhang           if (allrows[i]) rmap_i[j] = j;
1101ec3c739cSHong Zhang           else rmap_i[irow_i[j]] = j;
110217df9f7cSHong Zhang         }
110317df9f7cSHong Zhang #endif
1104059f6b74SHong Zhang       } else rmap[i] = NULL;
1105059f6b74SHong Zhang     }
110617df9f7cSHong Zhang 
110717df9f7cSHong Zhang     /* Update lens from offproc data */
110817df9f7cSHong Zhang     {
110917df9f7cSHong Zhang       PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i;
111017df9f7cSHong Zhang 
111117df9f7cSHong Zhang       ierr    = MPI_Waitall(nrqs,r_waits3,r_status3);CHKERRQ(ierr);
111217df9f7cSHong Zhang       for (tmp2=0; tmp2<nrqs; tmp2++) {
111317df9f7cSHong Zhang         sbuf1_i = sbuf1[pa[tmp2]];
111417df9f7cSHong Zhang         jmax    = sbuf1_i[0];
111517df9f7cSHong Zhang         ct1     = 2*jmax+1;
111617df9f7cSHong Zhang         ct2     = 0;
111717df9f7cSHong Zhang         rbuf2_i = rbuf2[tmp2];
111817df9f7cSHong Zhang         rbuf3_i = rbuf3[tmp2];
111917df9f7cSHong Zhang         for (j=1; j<=jmax; j++) {
112017df9f7cSHong Zhang           is_no  = sbuf1_i[2*j-1];
112117df9f7cSHong Zhang           max1   = sbuf1_i[2*j];
112217df9f7cSHong Zhang           lens_i = lens[is_no];
112317df9f7cSHong Zhang           if (!allcolumns[is_no]) cmap_i = cmap[is_no];
112417df9f7cSHong Zhang           rmap_i = rmap[is_no];
112517df9f7cSHong Zhang           for (k=0; k<max1; k++,ct1++) {
1126059f6b74SHong Zhang             if (allrows[is_no]) {
1127059f6b74SHong Zhang               row = sbuf1_i[ct1];
1128059f6b74SHong Zhang             } else {
112917df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
113017df9f7cSHong Zhang               ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
113117df9f7cSHong Zhang               row--;
113217df9f7cSHong Zhang               if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
113317df9f7cSHong Zhang #else
113417df9f7cSHong Zhang               row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
113517df9f7cSHong Zhang #endif
1136059f6b74SHong Zhang             }
113717df9f7cSHong Zhang             max2 = rbuf2_i[ct1];
113817df9f7cSHong Zhang             for (l=0; l<max2; l++,ct2++) {
113917df9f7cSHong Zhang               if (!allcolumns[is_no]) {
114017df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
114117df9f7cSHong Zhang                 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
114217df9f7cSHong Zhang #else
114317df9f7cSHong Zhang                 tcol = cmap_i[rbuf3_i[ct2]];
114417df9f7cSHong Zhang #endif
114517df9f7cSHong Zhang                 if (tcol) lens_i[row]++;
114617df9f7cSHong Zhang               } else { /* allcolumns */
114717df9f7cSHong Zhang                 lens_i[row]++; /* lens_i[row] += max2 ? */
114817df9f7cSHong Zhang               }
114917df9f7cSHong Zhang             }
115017df9f7cSHong Zhang           }
115117df9f7cSHong Zhang         }
115217df9f7cSHong Zhang       }
115317df9f7cSHong Zhang     }
115417df9f7cSHong Zhang     ierr = PetscFree(r_waits3);CHKERRQ(ierr);
115517df9f7cSHong Zhang     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
115617df9f7cSHong Zhang     ierr = PetscFree2(r_status3,s_status3);CHKERRQ(ierr);
115717df9f7cSHong Zhang     ierr = PetscFree(s_waits3);CHKERRQ(ierr);
115817df9f7cSHong Zhang 
115917df9f7cSHong Zhang     /* Create the submatrices */
116017df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
1161a42ab0d6SHong Zhang       PetscInt bs_tmp;
1162a42ab0d6SHong Zhang       if (ijonly) bs_tmp = 1;
1163a42ab0d6SHong Zhang       else        bs_tmp = bs;
116417df9f7cSHong Zhang 
116517df9f7cSHong Zhang       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
1166a42ab0d6SHong Zhang       ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
116717df9f7cSHong Zhang 
116817df9f7cSHong Zhang       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
1169a42ab0d6SHong Zhang       ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr);
1170a42ab0d6SHong Zhang       ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */
117117df9f7cSHong Zhang 
117217df9f7cSHong Zhang       /* create struct Mat_SubMat and attached it to submat */
117317df9f7cSHong Zhang       ierr = PetscNew(&smat_i);CHKERRQ(ierr);
117417df9f7cSHong Zhang       subc = (Mat_SeqBAIJ*)submats[i]->data;
117517df9f7cSHong Zhang       subc->submatis1 = smat_i;
117617df9f7cSHong Zhang       smats[i]        = smat_i;
117717df9f7cSHong Zhang 
117817df9f7cSHong Zhang       smat_i->destroy          = submats[i]->ops->destroy;
117917df9f7cSHong Zhang       submats[i]->ops->destroy = MatDestroy_MPIBAIJ_MatGetSubmatrices;
118017df9f7cSHong Zhang       submats[i]->factortype   = C->factortype;
118117df9f7cSHong Zhang 
118217df9f7cSHong Zhang       smat_i->id          = i;
118317df9f7cSHong Zhang       smat_i->nrqs        = nrqs;
118417df9f7cSHong Zhang       smat_i->nrqr        = nrqr;
118517df9f7cSHong Zhang       smat_i->rbuf1       = rbuf1;
118617df9f7cSHong Zhang       smat_i->rbuf2       = rbuf2;
118717df9f7cSHong Zhang       smat_i->rbuf3       = rbuf3;
118817df9f7cSHong Zhang       smat_i->sbuf2       = sbuf2;
118917df9f7cSHong Zhang       smat_i->req_source2 = req_source2;
119017df9f7cSHong Zhang 
119117df9f7cSHong Zhang       smat_i->sbuf1       = sbuf1;
119217df9f7cSHong Zhang       smat_i->ptr         = ptr;
119317df9f7cSHong Zhang       smat_i->tmp         = tmp;
119417df9f7cSHong Zhang       smat_i->ctr         = ctr;
119517df9f7cSHong Zhang 
119617df9f7cSHong Zhang       smat_i->pa           = pa;
119717df9f7cSHong Zhang       smat_i->req_size     = req_size;
119817df9f7cSHong Zhang       smat_i->req_source1  = req_source1;
119917df9f7cSHong Zhang 
120017df9f7cSHong Zhang       smat_i->allcolumns  = allcolumns[i];
1201ec3c739cSHong Zhang       smat_i->allrows     = allrows[i];
120217df9f7cSHong Zhang       smat_i->singleis    = PETSC_FALSE;
120317df9f7cSHong Zhang       smat_i->row2proc    = row2proc[i];
120417df9f7cSHong Zhang       smat_i->rmap        = rmap[i];
120517df9f7cSHong Zhang       smat_i->cmap        = cmap[i];
120617df9f7cSHong Zhang     }
120717df9f7cSHong Zhang 
120817df9f7cSHong Zhang     if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);}
120917df9f7cSHong Zhang     ierr = PetscFree(lens);CHKERRQ(ierr);
121017df9f7cSHong Zhang     ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
121117df9f7cSHong Zhang     ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
121217df9f7cSHong Zhang 
121317df9f7cSHong Zhang   } /* endof scall == MAT_INITIAL_MATRIX */
121417df9f7cSHong Zhang 
121517df9f7cSHong Zhang   /* Post recv matrix values */
1216bbc89d27SHong Zhang   if (!ijonly) {
121717df9f7cSHong Zhang     ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr);
121817df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr);
121917df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr);
122017df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr);
122117df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr);
122217df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
1223a42ab0d6SHong Zhang       ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr);
1224a42ab0d6SHong Zhang       ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRQ(ierr);
122517df9f7cSHong Zhang     }
122617df9f7cSHong Zhang 
122717df9f7cSHong Zhang     /* Allocate sending buffers for a->a, and send them off */
122817df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr);
122917df9f7cSHong Zhang     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
12309e686edcSHong Zhang 
1231a42ab0d6SHong Zhang     ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr);
1232a42ab0d6SHong Zhang     for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
123317df9f7cSHong Zhang 
123417df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr);
123517df9f7cSHong Zhang 
123617df9f7cSHong Zhang     for (i=0; i<nrqr; i++) {
123717df9f7cSHong Zhang       rbuf1_i   = rbuf1[i];
123817df9f7cSHong Zhang       sbuf_aa_i = sbuf_aa[i];
123917df9f7cSHong Zhang       ct1       = 2*rbuf1_i[0]+1;
124017df9f7cSHong Zhang       ct2       = 0;
124117df9f7cSHong Zhang       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
124217df9f7cSHong Zhang         kmax = rbuf1_i[2*j];
124317df9f7cSHong Zhang         for (k=0; k<kmax; k++,ct1++) {
124417df9f7cSHong Zhang           row    = rbuf1_i[ct1] - rstart;
1245a42ab0d6SHong Zhang           nzA    = a_i[row+1] - a_i[row];
1246a42ab0d6SHong Zhang           nzB    = b_i[row+1] - b_i[row];
124717df9f7cSHong Zhang           ncols  = nzA + nzB;
124817df9f7cSHong Zhang           cworkB = b_j + b_i[row];
1249a42ab0d6SHong Zhang           vworkA = a_a + a_i[row]*bs2;
1250a42ab0d6SHong Zhang           vworkB = b_a + b_i[row]*bs2;
125117df9f7cSHong Zhang 
125217df9f7cSHong Zhang           /* load the column values for this row into vals*/
1253a42ab0d6SHong Zhang           vals = sbuf_aa_i+ct2*bs2;
1254a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1255a42ab0d6SHong Zhang             if ((bmap[cworkB[l]]) < cstart) {
1256a42ab0d6SHong Zhang               ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1257a42ab0d6SHong Zhang             } else break;
1258a42ab0d6SHong Zhang           }
1259a42ab0d6SHong Zhang           imark = l;
1260a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1261bbc89d27SHong Zhang             ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1262a42ab0d6SHong Zhang           }
1263a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1264a42ab0d6SHong Zhang             ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1265a42ab0d6SHong Zhang           }
126617df9f7cSHong Zhang 
126717df9f7cSHong Zhang           ct2 += ncols;
126817df9f7cSHong Zhang         }
126917df9f7cSHong Zhang       }
1270a42ab0d6SHong Zhang       ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRQ(ierr);
127117df9f7cSHong Zhang     }
1272bbc89d27SHong Zhang   }
127317df9f7cSHong Zhang 
127417df9f7cSHong Zhang   if (!ismax) {
127517df9f7cSHong Zhang     ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr);
127617df9f7cSHong Zhang     ierr = PetscFree(rbuf1);CHKERRQ(ierr);
127717df9f7cSHong Zhang   }
127817df9f7cSHong Zhang 
127917df9f7cSHong Zhang   /* Assemble the matrices */
128017df9f7cSHong Zhang   /* First assemble the local rows */
128117df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
128217df9f7cSHong Zhang     row2proc_i = row2proc[i];
128317df9f7cSHong Zhang     subc      = (Mat_SeqBAIJ*)submats[i]->data;
128417df9f7cSHong Zhang     imat_ilen = subc->ilen;
128517df9f7cSHong Zhang     imat_j    = subc->j;
128617df9f7cSHong Zhang     imat_i    = subc->i;
128717df9f7cSHong Zhang     imat_a    = subc->a;
128817df9f7cSHong Zhang 
128917df9f7cSHong Zhang     if (!allcolumns[i]) cmap_i = cmap[i];
129017df9f7cSHong Zhang     rmap_i = rmap[i];
129117df9f7cSHong Zhang     irow_i = irow[i];
129217df9f7cSHong Zhang     jmax   = nrow[i];
129317df9f7cSHong Zhang     for (j=0; j<jmax; j++) {
1294ec3c739cSHong Zhang       if (allrows[i]) row = j;
1295ec3c739cSHong Zhang       else row  = irow_i[j];
129617df9f7cSHong Zhang       proc = row2proc_i[j];
1297a42ab0d6SHong Zhang 
129817df9f7cSHong Zhang       if (proc == rank) {
1299a42ab0d6SHong Zhang 
1300a42ab0d6SHong Zhang         row    = row - rstart;
1301a42ab0d6SHong Zhang         nzA    = a_i[row+1] - a_i[row];
1302a42ab0d6SHong Zhang         nzB    = b_i[row+1] - b_i[row];
1303a42ab0d6SHong Zhang         cworkA = a_j + a_i[row];
1304a42ab0d6SHong Zhang         cworkB = b_j + b_i[row];
1305a42ab0d6SHong Zhang         if (!ijonly) {
1306a42ab0d6SHong Zhang           vworkA = a_a + a_i[row]*bs2;
1307a42ab0d6SHong Zhang           vworkB = b_a + b_i[row]*bs2;
1308a42ab0d6SHong Zhang         }
1309059f6b74SHong Zhang 
1310059f6b74SHong Zhang         if (allrows[i]) {
1311059f6b74SHong Zhang           row = row+rstart;
1312059f6b74SHong Zhang         } else {
1313a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1314a42ab0d6SHong Zhang           ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr);
1315a42ab0d6SHong Zhang           row--;
1316121971b2SHong Zhang 
1317a42ab0d6SHong Zhang           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1318a42ab0d6SHong Zhang #else
1319a42ab0d6SHong Zhang           row = rmap_i[row + rstart];
1320a42ab0d6SHong Zhang #endif
1321059f6b74SHong Zhang         }
1322a42ab0d6SHong Zhang         mat_i = imat_i[row];
1323a42ab0d6SHong Zhang         if (!ijonly) mat_a = imat_a + mat_i*bs2;
1324a42ab0d6SHong Zhang         mat_j    = imat_j + mat_i;
132580d31651SHong Zhang         ilen = imat_ilen[row];
1326a42ab0d6SHong Zhang 
1327a42ab0d6SHong Zhang         /* load the column indices for this row into cols*/
1328a42ab0d6SHong Zhang         if (!allcolumns[i]) {
1329a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1330a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1331a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1332a42ab0d6SHong Zhang               ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr);
1333a42ab0d6SHong Zhang               if (tcol) {
1334a42ab0d6SHong Zhang #else
1335a42ab0d6SHong Zhang               if ((tcol = cmap_i[ctmp])) {
1336a42ab0d6SHong Zhang #endif
1337a42ab0d6SHong Zhang                 *mat_j++ = tcol - 1;
1338a42ab0d6SHong Zhang                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1339a42ab0d6SHong Zhang                 mat_a   += bs2;
134080d31651SHong Zhang                 ilen++;
1341a42ab0d6SHong Zhang               }
1342a42ab0d6SHong Zhang             } else break;
1343a42ab0d6SHong Zhang           }
1344a42ab0d6SHong Zhang           imark = l;
1345a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1346a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1347a42ab0d6SHong Zhang             ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1348a42ab0d6SHong Zhang             if (tcol) {
1349a42ab0d6SHong Zhang #else
1350a42ab0d6SHong Zhang             if ((tcol = cmap_i[cstart + cworkA[l]])) {
1351a42ab0d6SHong Zhang #endif
1352a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1353a42ab0d6SHong Zhang               if (!ijonly) {
1354a42ab0d6SHong Zhang                 ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1355a42ab0d6SHong Zhang                 mat_a += bs2;
1356a42ab0d6SHong Zhang               }
135780d31651SHong Zhang               ilen++;
1358a42ab0d6SHong Zhang             }
1359a42ab0d6SHong Zhang           }
1360a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1361a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1362a42ab0d6SHong Zhang             ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1363a42ab0d6SHong Zhang             if (tcol) {
1364a42ab0d6SHong Zhang #else
1365a42ab0d6SHong Zhang             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1366a42ab0d6SHong Zhang #endif
1367a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1368a42ab0d6SHong Zhang               if (!ijonly) {
1369a42ab0d6SHong Zhang                 ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1370a42ab0d6SHong Zhang                 mat_a += bs2;
1371a42ab0d6SHong Zhang               }
137280d31651SHong Zhang               ilen++;
1373a42ab0d6SHong Zhang             }
1374a42ab0d6SHong Zhang           }
1375a42ab0d6SHong Zhang         } else { /* allcolumns */
1376a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1377a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1378a42ab0d6SHong Zhang               *mat_j++ = ctmp;
1379a42ab0d6SHong Zhang               ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1380a42ab0d6SHong Zhang               mat_a   += bs2;
138180d31651SHong Zhang               ilen++;
1382a42ab0d6SHong Zhang             } else break;
1383a42ab0d6SHong Zhang           }
1384a42ab0d6SHong Zhang           imark = l;
1385a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1386a42ab0d6SHong Zhang             *mat_j++ = cstart+cworkA[l];
1387a42ab0d6SHong Zhang             if (!ijonly) {
1388a42ab0d6SHong Zhang               ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1389a42ab0d6SHong Zhang               mat_a += bs2;
1390a42ab0d6SHong Zhang             }
139180d31651SHong Zhang             ilen++;
1392a42ab0d6SHong Zhang           }
1393a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1394a42ab0d6SHong Zhang             *mat_j++ = bmap[cworkB[l]];
1395a42ab0d6SHong Zhang             if (!ijonly) {
1396a42ab0d6SHong Zhang               ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1397a42ab0d6SHong Zhang               mat_a += bs2;
1398a42ab0d6SHong Zhang             }
139980d31651SHong Zhang             ilen++;
1400a42ab0d6SHong Zhang           }
1401a42ab0d6SHong Zhang         }
140280d31651SHong Zhang         imat_ilen[row] = ilen;
140317df9f7cSHong Zhang       }
140417df9f7cSHong Zhang     }
140517df9f7cSHong Zhang   }
140617df9f7cSHong Zhang 
140717df9f7cSHong Zhang   /* Now assemble the off proc rows */
14085dd0c08cSHong Zhang   if (!ijonly) {
140917df9f7cSHong Zhang     ierr = MPI_Waitall(nrqs,r_waits4,r_status4);CHKERRQ(ierr);
14105dd0c08cSHong Zhang   }
141117df9f7cSHong Zhang   for (tmp2=0; tmp2<nrqs; tmp2++) {
141217df9f7cSHong Zhang     sbuf1_i = sbuf1[pa[tmp2]];
141317df9f7cSHong Zhang     jmax    = sbuf1_i[0];
141417df9f7cSHong Zhang     ct1     = 2*jmax + 1;
141517df9f7cSHong Zhang     ct2     = 0;
141617df9f7cSHong Zhang     rbuf2_i = rbuf2[tmp2];
141717df9f7cSHong Zhang     rbuf3_i = rbuf3[tmp2];
14185dd0c08cSHong Zhang     if (!ijonly) rbuf4_i = rbuf4[tmp2];
141917df9f7cSHong Zhang     for (j=1; j<=jmax; j++) {
142017df9f7cSHong Zhang       is_no     = sbuf1_i[2*j-1];
142117df9f7cSHong Zhang       rmap_i    = rmap[is_no];
142217df9f7cSHong Zhang       if (!allcolumns[is_no]) cmap_i = cmap[is_no];
142317df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ*)submats[is_no]->data;
142417df9f7cSHong Zhang       imat_ilen = subc->ilen;
142517df9f7cSHong Zhang       imat_j    = subc->j;
142617df9f7cSHong Zhang       imat_i    = subc->i;
14275dd0c08cSHong Zhang       if (!ijonly) imat_a    = subc->a;
142817df9f7cSHong Zhang       max1      = sbuf1_i[2*j];
14295dd0c08cSHong Zhang       for (k=0; k<max1; k++,ct1++) { /* for each recved block row */
143017df9f7cSHong Zhang         row = sbuf1_i[ct1];
1431059f6b74SHong Zhang 
1432059f6b74SHong Zhang         if (allrows[is_no]) {
1433059f6b74SHong Zhang           row = sbuf1_i[ct1];
1434059f6b74SHong Zhang         } else {
143517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
143617df9f7cSHong Zhang           ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr);
143717df9f7cSHong Zhang           row--;
14385dd0c08cSHong Zhang           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
143917df9f7cSHong Zhang #else
144017df9f7cSHong Zhang           row = rmap_i[row];
144117df9f7cSHong Zhang #endif
1442059f6b74SHong Zhang         }
144317df9f7cSHong Zhang         ilen  = imat_ilen[row];
144417df9f7cSHong Zhang         mat_i = imat_i[row];
14455dd0c08cSHong Zhang         if (!ijonly) mat_a = imat_a + mat_i*bs2;
144617df9f7cSHong Zhang         mat_j = imat_j + mat_i;
144717df9f7cSHong Zhang         max2  = rbuf2_i[ct1];
144817df9f7cSHong Zhang         if (!allcolumns[is_no]) {
144917df9f7cSHong Zhang           for (l=0; l<max2; l++,ct2++) {
145017df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
145117df9f7cSHong Zhang             ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
145217df9f7cSHong Zhang #else
145317df9f7cSHong Zhang             tcol = cmap_i[rbuf3_i[ct2]];
145417df9f7cSHong Zhang #endif
145517df9f7cSHong Zhang             if (tcol) {
145617df9f7cSHong Zhang               *mat_j++ = tcol - 1;
14575dd0c08cSHong Zhang               if (!ijonly) {
14585dd0c08cSHong Zhang                 ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
14595dd0c08cSHong Zhang                 mat_a += bs2;
14605dd0c08cSHong Zhang               }
146117df9f7cSHong Zhang               ilen++;
146217df9f7cSHong Zhang             }
146317df9f7cSHong Zhang           }
146417df9f7cSHong Zhang         } else { /* allcolumns */
146517df9f7cSHong Zhang           for (l=0; l<max2; l++,ct2++) {
146617df9f7cSHong Zhang             *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */
14675dd0c08cSHong Zhang             if (!ijonly) {
14685dd0c08cSHong Zhang               ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
14695dd0c08cSHong Zhang               mat_a += bs2;
14705dd0c08cSHong Zhang             }
147117df9f7cSHong Zhang             ilen++;
147217df9f7cSHong Zhang           }
147317df9f7cSHong Zhang         }
147417df9f7cSHong Zhang         imat_ilen[row] = ilen;
147517df9f7cSHong Zhang       }
147617df9f7cSHong Zhang     }
147717df9f7cSHong Zhang   }
147817df9f7cSHong Zhang 
147980d31651SHong Zhang   if (!iscsorted) { /* sort column indices of the rows */
148080d31651SHong Zhang     MatScalar *work;
148180d31651SHong Zhang 
148280d31651SHong Zhang     ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr);
148317df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
148417df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ*)submats[i]->data;
148580d31651SHong Zhang       imat_ilen = subc->ilen;
148617df9f7cSHong Zhang       imat_j    = subc->j;
148717df9f7cSHong Zhang       imat_i    = subc->i;
14884b6ceb0dSHong Zhang       if (!ijonly) imat_a = subc->a;
148917df9f7cSHong Zhang       if (allcolumns[i]) continue;
14904b6ceb0dSHong Zhang 
149117df9f7cSHong Zhang       jmax = nrow[i];
149217df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
149380d31651SHong Zhang         mat_i = imat_i[j];
149480d31651SHong Zhang         mat_j = imat_j + mat_i;
14954b6ceb0dSHong Zhang         ilen  = imat_ilen[j];
149680d31651SHong Zhang         if (ijonly) {
149780d31651SHong Zhang           ierr = PetscSortInt(ilen,mat_j);CHKERRQ(ierr);
149880d31651SHong Zhang         } else {
14994b6ceb0dSHong Zhang           mat_a = imat_a + mat_i*bs2;
150080d31651SHong Zhang           ierr = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);
150117df9f7cSHong Zhang         }
150217df9f7cSHong Zhang       }
150317df9f7cSHong Zhang     }
150480d31651SHong Zhang     ierr = PetscFree(work);CHKERRQ(ierr);
150580d31651SHong Zhang   }
150617df9f7cSHong Zhang 
1507bbc89d27SHong Zhang   if (!ijonly) {
150817df9f7cSHong Zhang     ierr = PetscFree(r_status4);CHKERRQ(ierr);
150917df9f7cSHong Zhang     ierr = PetscFree(r_waits4);CHKERRQ(ierr);
151017df9f7cSHong Zhang     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
151117df9f7cSHong Zhang     ierr = PetscFree(s_waits4);CHKERRQ(ierr);
151217df9f7cSHong Zhang     ierr = PetscFree(s_status4);CHKERRQ(ierr);
1513bbc89d27SHong Zhang   }
151417df9f7cSHong Zhang 
151517df9f7cSHong Zhang   /* Restore the indices */
151617df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
1517ec3c739cSHong Zhang     if (!allrows[i]) {
151817df9f7cSHong Zhang       ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1519ec3c739cSHong Zhang     }
152017df9f7cSHong Zhang     if (!allcolumns[i]) {
152117df9f7cSHong Zhang       ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
152217df9f7cSHong Zhang     }
152317df9f7cSHong Zhang   }
152417df9f7cSHong Zhang 
152517df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
152617df9f7cSHong Zhang     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
152717df9f7cSHong Zhang     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
152817df9f7cSHong Zhang   }
152917df9f7cSHong Zhang 
153017df9f7cSHong Zhang   /* Destroy allocated memory */
153117df9f7cSHong Zhang   if (!ismax) {
153217df9f7cSHong Zhang     ierr = PetscFree(pa);CHKERRQ(ierr);
153317df9f7cSHong Zhang     ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
153417df9f7cSHong Zhang     for (i=0; i<nrqr; ++i) {
153517df9f7cSHong Zhang       ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
153617df9f7cSHong Zhang     }
153717df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
153817df9f7cSHong Zhang       ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
153917df9f7cSHong Zhang     }
154017df9f7cSHong Zhang 
154117df9f7cSHong Zhang     ierr = PetscFree3(sbuf2,req_size,req_source1);CHKERRQ(ierr);
154217df9f7cSHong Zhang     ierr = PetscFree3(req_source2,rbuf2,rbuf3);CHKERRQ(ierr);
154317df9f7cSHong Zhang   }
1544bbc89d27SHong Zhang   ierr = PetscFree5(irow,icol,nrow,ncol,issorted);CHKERRQ(ierr);
1545bbc89d27SHong Zhang   ierr = PetscFree6(smats,row2proc,cmap,rmap,allcolumns,allrows);CHKERRQ(ierr);
154617df9f7cSHong Zhang 
1547bbc89d27SHong Zhang   if (!ijonly) {
154817df9f7cSHong Zhang     ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
154917df9f7cSHong Zhang     ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
155017df9f7cSHong Zhang 
155117df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
155217df9f7cSHong Zhang       ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);
155317df9f7cSHong Zhang     }
155417df9f7cSHong Zhang     ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1555bbc89d27SHong Zhang   }
15567a868f3eSHong Zhang   c->ijonly = PETSC_FALSE; /* set back to the default */
15573a40ed3dSBarry Smith   PetscFunctionReturn(0);
1558d5b485abSSatish Balay }
1559