xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 4e195584755f7b76b0f8e73df02b8ccefab61707)
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)
562d4ee042Sprj-   nrqr - no of requests received (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;
65131c27b5Sprj-   PetscInt       Mbs,i,j,k,**rbuf,row,nrqs,msz,**outdat,**ptr;
668f9f447aSHong Zhang   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p;
67131c27b5Sprj-   PetscMPIInt    *onodes1,*olengths1,*onodes2,*olengths2,proc=-1;
68f1af5d2fSBarry Smith   PetscBT        *table;
69749130a8SStefano Zampini   MPI_Comm       comm,*iscomms;
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 
83f489ac74SBarry Smith   ierr = PetscMalloc2(imax+1,(PetscInt***)&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++) {
94580bdb30SBarry Smith     ierr  = PetscArrayzero(w4,size);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);
137580bdb30SBarry Smith   ierr = PetscArrayzero(outdat,size);CHKERRQ(ierr);
138580bdb30SBarry Smith   ierr = PetscArrayzero(ptr,size);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;
154580bdb30SBarry Smith     ierr         = PetscArrayzero(outdat[j]+1,2*w3[j]);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++) {
174580bdb30SBarry Smith       ierr    = PetscArrayzero(ctr,size);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];
208ffc4695bSBarry Smith     ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRMPI(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   }
2152cd485c2SBarry Smith   ierr = PetscFree2(*(PetscInt***)&idx,n);CHKERRQ(ierr);
216d5b485abSSatish Balay 
217749130a8SStefano Zampini   ierr = PetscMalloc1(imax,&iscomms);CHKERRQ(ierr);
218d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
219749130a8SStefano Zampini     ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]),&iscomms[i],NULL);CHKERRQ(ierr);
2206bf464f9SBarry Smith     ierr = ISDestroy(&is[i]);CHKERRQ(ierr);
221d5b485abSSatish Balay   }
222d5b485abSSatish Balay 
223d5b485abSSatish Balay   /* Do Local work*/
224df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
225d5b485abSSatish Balay 
226d5b485abSSatish Balay   /* Receive messages*/
227854ce69bSBarry Smith   ierr = PetscMalloc1(nrqr+1,&recv_status);CHKERRQ(ierr);
228ffc4695bSBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRMPI(ierr);}
229d5b485abSSatish Balay 
230854ce69bSBarry Smith   ierr = PetscMalloc1(nrqs+1,&s_status);CHKERRQ(ierr);
231ffc4695bSBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRMPI(ierr);}
232d5b485abSSatish Balay 
233d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
23423bdbc58SBarry Smith   ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr);
23523bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
236d5b485abSSatish Balay 
237854ce69bSBarry Smith   ierr = PetscMalloc1(nrqr+1,&xdata);CHKERRQ(ierr);
238854ce69bSBarry Smith   ierr = PetscMalloc1(nrqr+1,&isz1);CHKERRQ(ierr);
239df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
240c05d87d6SBarry Smith   ierr = PetscFree(rbuf[0]);CHKERRQ(ierr);
241606d414cSSatish Balay   ierr = PetscFree(rbuf);CHKERRQ(ierr);
242d5b485abSSatish Balay 
243d5b485abSSatish Balay   /* Send the data back*/
244d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
245d5b485abSSatish Balay   {
246b24ad042SBarry Smith     PetscMPIInt *rw1;
247d5b485abSSatish Balay 
248884e879aSBarry Smith     ierr = PetscCalloc1(size,&rw1);CHKERRQ(ierr);
249c7dd2924SSatish Balay 
250d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
251d5b485abSSatish Balay       proc = recv_status[i].MPI_SOURCE;
252e32f2f54SBarry Smith       if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch");
253d5b485abSSatish Balay       rw1[proc] = isz1[i];
254d5b485abSSatish Balay     }
255d5b485abSSatish Balay 
256c7dd2924SSatish Balay     ierr = PetscFree(onodes1);CHKERRQ(ierr);
257c7dd2924SSatish Balay     ierr = PetscFree(olengths1);CHKERRQ(ierr);
258d5b485abSSatish Balay 
259c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
260c7dd2924SSatish Balay     ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr);
26103512dc5SSatish Balay     ierr = PetscFree(rw1);CHKERRQ(ierr);
262c7dd2924SSatish Balay   }
263c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
264c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr);
265d5b485abSSatish Balay 
266d5b485abSSatish Balay   /*  Now  post the sends */
267854ce69bSBarry Smith   ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr);
268d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
269d5b485abSSatish Balay     j    = recv_status[i].MPI_SOURCE;
270ffc4695bSBarry Smith     ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRMPI(ierr);
271d5b485abSSatish Balay   }
272d5b485abSSatish Balay 
273d5b485abSSatish Balay   /* receive work done on other processors*/
274d5b485abSSatish Balay   {
275b24ad042SBarry Smith     PetscMPIInt idex;
276b24ad042SBarry Smith     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
277f1af5d2fSBarry Smith     PetscBT     table_i;
278d5b485abSSatish Balay     MPI_Status  *status2;
279d5b485abSSatish Balay 
280854ce69bSBarry Smith     ierr = PetscMalloc1(PetscMax(nrqr,nrqs)+1,&status2);CHKERRQ(ierr);
281d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
282ffc4695bSBarry Smith       ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRMPI(ierr);
283d5b485abSSatish Balay       /* Process the message*/
284999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
285d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
286999d9058SBarry Smith       jmax    = rbuf2[idex][0];
287d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
288d5b485abSSatish Balay         max     = rbuf2_i[2*j];
289d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
290d5b485abSSatish Balay         isz_i   = isz[is_no];
291d5b485abSSatish Balay         data_i  = data[is_no];
292d5b485abSSatish Balay         table_i = table[is_no];
293d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
294d5b485abSSatish Balay           row = rbuf2_i[ct1];
29526fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row;
296d5b485abSSatish Balay         }
297d5b485abSSatish Balay         isz[is_no] = isz_i;
298d5b485abSSatish Balay       }
299d5b485abSSatish Balay     }
300ffc4695bSBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRMPI(ierr);}
301606d414cSSatish Balay     ierr = PetscFree(status2);CHKERRQ(ierr);
302d5b485abSSatish Balay   }
303d5b485abSSatish Balay 
304d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
305749130a8SStefano Zampini     ierr = ISCreateGeneral(iscomms[i],isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
306749130a8SStefano Zampini     ierr = PetscCommDestroy(&iscomms[i]);CHKERRQ(ierr);
307d5b485abSSatish Balay   }
308d5b485abSSatish Balay 
309749130a8SStefano Zampini   ierr = PetscFree(iscomms);CHKERRQ(ierr);
310c7dd2924SSatish Balay   ierr = PetscFree(onodes2);CHKERRQ(ierr);
311c7dd2924SSatish Balay   ierr = PetscFree(olengths2);CHKERRQ(ierr);
312c7dd2924SSatish Balay 
313606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
314c05d87d6SBarry Smith   ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr);
315606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
316606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
317606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
318606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
319606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
3208f9f447aSHong Zhang   ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr);
321606d414cSSatish Balay   ierr = PetscFree(s_status);CHKERRQ(ierr);
322606d414cSSatish Balay   ierr = PetscFree(recv_status);CHKERRQ(ierr);
323606d414cSSatish Balay   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
324606d414cSSatish Balay   ierr = PetscFree(xdata);CHKERRQ(ierr);
325606d414cSSatish Balay   ierr = PetscFree(isz1);CHKERRQ(ierr);
3263a40ed3dSBarry Smith   PetscFunctionReturn(0);
327d5b485abSSatish Balay }
328d5b485abSSatish Balay 
329d5b485abSSatish Balay /*
330df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
331d5b485abSSatish Balay        the work on the local processor.
332d5b485abSSatish Balay 
333d5b485abSSatish Balay      Inputs:
334df36731bSBarry Smith       C      - MAT_MPIBAIJ;
335d5b485abSSatish Balay       imax - total no of index sets processed at a time;
336df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
337d5b485abSSatish Balay 
338d5b485abSSatish Balay      Output:
3390e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
340d5b485abSSatish Balay                to each index set;
341d5b485abSSatish Balay       data   - pointer to the solutions
342d5b485abSSatish Balay */
343b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data)
344d5b485abSSatish Balay {
345df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
346d5b485abSSatish Balay   Mat         A  = c->A,B = c->B;
347df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
348b24ad042SBarry Smith   PetscInt    start,end,val,max,rstart,cstart,*ai,*aj;
349b24ad042SBarry Smith   PetscInt    *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
350f1af5d2fSBarry Smith   PetscBT     table_i;
351d5b485abSSatish Balay 
3523a40ed3dSBarry Smith   PetscFunctionBegin;
353899cda47SBarry Smith   rstart = c->rstartbs;
354899cda47SBarry Smith   cstart = c->cstartbs;
355d5b485abSSatish Balay   ai     = a->i;
356df36731bSBarry Smith   aj     = a->j;
357d5b485abSSatish Balay   bi     = b->i;
358df36731bSBarry Smith   bj     = b->j;
359d5b485abSSatish Balay   garray = c->garray;
360d5b485abSSatish Balay 
361d5b485abSSatish Balay 
362d5b485abSSatish Balay   for (i=0; i<imax; i++) {
363d5b485abSSatish Balay     data_i  = data[i];
364d5b485abSSatish Balay     table_i = table[i];
365d5b485abSSatish Balay     isz_i   = isz[i];
366d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
367d5b485abSSatish Balay       row   = data_i[j] - rstart;
368d5b485abSSatish Balay       start = ai[row];
369d5b485abSSatish Balay       end   = ai[row+1];
370d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
371df36731bSBarry Smith         val = aj[k] + cstart;
37226fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val;
373d5b485abSSatish Balay       }
374d5b485abSSatish Balay       start = bi[row];
375d5b485abSSatish Balay       end   = bi[row+1];
376d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
377df36731bSBarry Smith         val = garray[bj[k]];
37826fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val;
379d5b485abSSatish Balay       }
380d5b485abSSatish Balay     }
381d5b485abSSatish Balay     isz[i] = isz_i;
382d5b485abSSatish Balay   }
3833a40ed3dSBarry Smith   PetscFunctionReturn(0);
384d5b485abSSatish Balay }
385d5b485abSSatish Balay /*
3862d4ee042Sprj-       MatIncreaseOverlap_MPIBAIJ_Receive - Process the received messages,
387d5b485abSSatish Balay          and return the output
388d5b485abSSatish Balay 
389d5b485abSSatish Balay          Input:
390d5b485abSSatish Balay            C    - the matrix
391d5b485abSSatish Balay            nrqr - no of messages being processed.
3922d4ee042Sprj-            rbuf - an array of pointers to the received requests
393d5b485abSSatish Balay 
394d5b485abSSatish Balay          Output:
395d5b485abSSatish Balay            xdata - array of messages to be sent back
396d5b485abSSatish Balay            isz1  - size of each message
397d5b485abSSatish Balay 
398a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
399d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
4000e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of
401d5b485abSSatish Balay memory is used.
402d5b485abSSatish Balay 
403d5b485abSSatish Balay */
404b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
405d5b485abSSatish Balay {
406df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
407d5b485abSSatish Balay   Mat            A  = c->A,B = c->B;
408df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
4096849ba73SBarry Smith   PetscErrorCode ierr;
410b24ad042SBarry Smith   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
411b24ad042SBarry Smith   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
412d2a212eaSBarry Smith   PetscInt       val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
413b24ad042SBarry Smith   PetscInt       *rbuf_i,kmax,rbuf_0;
414f1af5d2fSBarry Smith   PetscBT        xtable;
415d5b485abSSatish Balay 
4163a40ed3dSBarry Smith   PetscFunctionBegin;
417df36731bSBarry Smith   Mbs    = c->Mbs;
418899cda47SBarry Smith   rstart = c->rstartbs;
419899cda47SBarry Smith   cstart = c->cstartbs;
420d5b485abSSatish Balay   ai     = a->i;
421df36731bSBarry Smith   aj     = a->j;
422d5b485abSSatish Balay   bi     = b->i;
423df36731bSBarry Smith   bj     = b->j;
424d5b485abSSatish Balay   garray = c->garray;
425d5b485abSSatish Balay 
426d5b485abSSatish Balay 
427d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
428d5b485abSSatish Balay     rbuf_i =  rbuf[i];
429d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
430d5b485abSSatish Balay     ct    += rbuf_0;
43126fbe8dcSKarl Rupp     for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j];
432d5b485abSSatish Balay   }
433d5b485abSSatish Balay 
434701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
435701b8814SBarry Smith   else        max1 = 1;
436d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
437785e854fSJed Brown   ierr         = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr);
438d5b485abSSatish Balay   ++no_malloc;
43953b8de81SBarry Smith   ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr);
440580bdb30SBarry Smith   ierr = PetscArrayzero(isz1,nrqr);CHKERRQ(ierr);
441d5b485abSSatish Balay 
442d5b485abSSatish Balay   ct3 = 0;
443d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
444d5b485abSSatish Balay     rbuf_i =  rbuf[i];
445d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
446d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
447d5b485abSSatish Balay     ct2    =  ct1;
448d5b485abSSatish Balay     ct3   += ct1;
449d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
4506831982aSBarry Smith       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
451d5b485abSSatish Balay       oct2 = ct2;
452d5b485abSSatish Balay       kmax = rbuf_i[2*j];
453d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
454d5b485abSSatish Balay         row = rbuf_i[ct1];
455f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
456d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
457b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
458785e854fSJed Brown             ierr         = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr);
459580bdb30SBarry Smith             ierr         = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr);
460606d414cSSatish Balay             ierr         = PetscFree(xdata[0]);CHKERRQ(ierr);
461d5b485abSSatish Balay             xdata[0]     = tmp;
462d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
46326fbe8dcSKarl Rupp             for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
464d5b485abSSatish Balay           }
465d5b485abSSatish Balay           xdata[i][ct2++] = row;
466d5b485abSSatish Balay           ct3++;
467d5b485abSSatish Balay         }
468d5b485abSSatish Balay       }
469d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
470d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
471d5b485abSSatish Balay         start = ai[row];
472d5b485abSSatish Balay         end   = ai[row+1];
473d5b485abSSatish Balay         for (l=start; l<end; l++) {
474df36731bSBarry Smith           val = aj[l] + cstart;
475f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
476d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
477b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
478785e854fSJed Brown               ierr         = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr);
479580bdb30SBarry Smith               ierr         = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr);
480606d414cSSatish Balay               ierr         = PetscFree(xdata[0]);CHKERRQ(ierr);
481d5b485abSSatish Balay               xdata[0]     = tmp;
482d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
48326fbe8dcSKarl Rupp               for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
484d5b485abSSatish Balay             }
485d5b485abSSatish Balay             xdata[i][ct2++] = val;
486d5b485abSSatish Balay             ct3++;
487d5b485abSSatish Balay           }
488d5b485abSSatish Balay         }
489d5b485abSSatish Balay         start = bi[row];
490d5b485abSSatish Balay         end   = bi[row+1];
491d5b485abSSatish Balay         for (l=start; l<end; l++) {
492df36731bSBarry Smith           val = garray[bj[l]];
493f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
494d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
495b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
496785e854fSJed Brown               ierr         = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr);
497580bdb30SBarry Smith               ierr         = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr);
498606d414cSSatish Balay               ierr         = PetscFree(xdata[0]);CHKERRQ(ierr);
499d5b485abSSatish Balay               xdata[0]     = tmp;
500d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
50126fbe8dcSKarl Rupp               for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
502d5b485abSSatish Balay             }
503d5b485abSSatish Balay             xdata[i][ct2++] = val;
504d5b485abSSatish Balay             ct3++;
505d5b485abSSatish Balay           }
506d5b485abSSatish Balay         }
507d5b485abSSatish Balay       }
508d5b485abSSatish Balay       /* Update the header*/
509d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
510d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
511d5b485abSSatish Balay     }
512d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
513d5b485abSSatish Balay     xdata[i+1]  = xdata[i] + ct2;
514d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
515d5b485abSSatish Balay   }
51694bacf5dSBarry Smith   ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr);
5171e2582c4SBarry Smith   ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr);
5183a40ed3dSBarry Smith   PetscFunctionReturn(0);
519d5b485abSSatish Balay }
520d5b485abSSatish Balay 
5217dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
522d5b485abSSatish Balay {
52317df9f7cSHong Zhang   IS             *isrow_block,*iscol_block;
524cf4f527aSSatish Balay   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
5256849ba73SBarry Smith   PetscErrorCode ierr;
526121971b2SHong Zhang   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,N=C->cmap->N,bs=C->rmap->bs;
5275dc98d6aSHong Zhang   Mat_SeqBAIJ    *subc;
5285c39f6d9SHong Zhang   Mat_SubSppt    *smat;
529a2feddc5SSatish Balay 
5303a40ed3dSBarry Smith   PetscFunctionBegin;
53129dcf524SDmitry Karpeev   /* The compression and expansion should be avoided. Doesn't point
53229dcf524SDmitry Karpeev      out errors, might change the indices, hence buggey */
53317df9f7cSHong Zhang   ierr = PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block);CHKERRQ(ierr);
53417df9f7cSHong Zhang   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block);CHKERRQ(ierr);
53517df9f7cSHong Zhang   ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block);CHKERRQ(ierr);
536cf4f527aSSatish Balay 
537cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
5384698041cSHong Zhang   if (!C->cmap->N) nmax=20*1000000/sizeof(PetscInt);
5394698041cSHong Zhang   else nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt));
540329f5518SBarry Smith   if (!nmax) nmax = 1;
5415dc98d6aSHong Zhang 
5425dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
5434698041cSHong Zhang     nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */
544cf4f527aSSatish Balay 
545653e4784SBarry Smith     /* Make sure every processor loops through the nstages */
546b2566f29SBarry Smith     ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr);
5475dc98d6aSHong Zhang 
5484698041cSHong Zhang     /* Allocate memory to hold all the submatrices and dummy submatrices */
5494698041cSHong Zhang     ierr = PetscCalloc1(ismax+nstages,submat);CHKERRQ(ierr);
5504698041cSHong Zhang   } else { /* MAT_REUSE_MATRIX */
5514698041cSHong Zhang     if (ismax) {
5525dc98d6aSHong Zhang       subc = (Mat_SeqBAIJ*)((*submat)[0]->data);
5535dc98d6aSHong Zhang       smat   = subc->submatis1;
5544698041cSHong Zhang     } else { /* (*submat)[0] is a dummy matrix */
5555c39f6d9SHong Zhang       smat = (Mat_SubSppt*)(*submat)[0]->data;
5564698041cSHong Zhang     }
557071fcb05SBarry Smith     if (!smat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat");
5585dc98d6aSHong Zhang     nstages = smat->nstages;
5595dc98d6aSHong Zhang   }
5605dc98d6aSHong Zhang 
561cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
562cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
563*4e195584SJunchao Zhang     else if (pos >= ismax) max_no = 0;
564cf4f527aSSatish Balay     else                   max_no = ismax-pos;
565a42ab0d6SHong Zhang 
5667dae84e0SHong Zhang     ierr = MatCreateSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,*submat+pos);CHKERRQ(ierr);
567*4e195584SJunchao Zhang     if (!max_no) {
568*4e195584SJunchao Zhang       if (scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */
5695c39f6d9SHong Zhang         smat = (Mat_SubSppt*)(*submat)[pos]->data;
5704698041cSHong Zhang         smat->nstages = nstages;
5714698041cSHong Zhang       }
572*4e195584SJunchao Zhang       pos++; /* advance to next dummy matrix if any */
573*4e195584SJunchao Zhang     } else pos += max_no;
574cf4f527aSSatish Balay   }
57536f4e84dSSatish Balay 
5765dc98d6aSHong Zhang   if (scall == MAT_INITIAL_MATRIX && ismax) {
5775dc98d6aSHong Zhang     /* save nstages for reuse */
5785dc98d6aSHong Zhang     subc = (Mat_SeqBAIJ*)((*submat)[0]->data);
5795dc98d6aSHong Zhang     smat = subc->submatis1;
5805dc98d6aSHong Zhang     smat->nstages = nstages;
5815dc98d6aSHong Zhang   }
5825dc98d6aSHong Zhang 
58336f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
58417df9f7cSHong Zhang     ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr);
58517df9f7cSHong Zhang     ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr);
58636f4e84dSSatish Balay   }
58717df9f7cSHong Zhang   ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr);
5883a40ed3dSBarry Smith   PetscFunctionReturn(0);
589a2feddc5SSatish Balay }
590a2feddc5SSatish Balay 
591233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
592e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
593233c2ff6SSatish Balay {
594e005ede5SBarry Smith   PetscInt       nGlobalNd = proc_gnode[size];
5954dc2109aSBarry Smith   PetscMPIInt    fproc;
5964dc2109aSBarry Smith   PetscErrorCode ierr;
597233c2ff6SSatish Balay 
598233c2ff6SSatish Balay   PetscFunctionBegin;
5994dc2109aSBarry Smith   ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr);
60023ce1328SBarry Smith   if (fproc > size) fproc = size;
601e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
602e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
603233c2ff6SSatish Balay     else                         fproc++;
604233c2ff6SSatish Balay   }
605e005ede5SBarry Smith   *rank = fproc;
606233c2ff6SSatish Balay   PetscFunctionReturn(0);
607233c2ff6SSatish Balay }
608233c2ff6SSatish Balay #endif
609233c2ff6SSatish Balay 
610a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
611b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
6127dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
61317df9f7cSHong Zhang {
61417df9f7cSHong Zhang   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
61517df9f7cSHong Zhang   Mat            A  = c->A;
61617df9f7cSHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc;
61717df9f7cSHong Zhang   const PetscInt **icol,**irow;
61817df9f7cSHong Zhang   PetscInt       *nrow,*ncol,start;
61917df9f7cSHong Zhang   PetscErrorCode ierr;
62017df9f7cSHong Zhang   PetscMPIInt    rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr;
621ddb3d6bcSHong Zhang   PetscInt       **sbuf1,**sbuf2,*sbuf2_i,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1;
62217df9f7cSHong Zhang   PetscInt       nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol;
62317df9f7cSHong Zhang   PetscInt       **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2;
62417df9f7cSHong Zhang   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
62517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
62617df9f7cSHong Zhang   PetscTable     *cmap,cmap_i=NULL,*rmap,rmap_i;
62717df9f7cSHong Zhang #else
62817df9f7cSHong Zhang   PetscInt       **cmap,*cmap_i=NULL,**rmap,*rmap_i;
62917df9f7cSHong Zhang #endif
63096cff833SHong Zhang   const PetscInt *irow_i,*icol_i;
63117df9f7cSHong Zhang   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
63217df9f7cSHong Zhang   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
63317df9f7cSHong Zhang   MPI_Request    *r_waits4,*s_waits3,*s_waits4;
63417df9f7cSHong Zhang   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
63517df9f7cSHong Zhang   MPI_Status     *r_status3,*r_status4,*s_status4;
63617df9f7cSHong Zhang   MPI_Comm       comm;
63764f5bb2dSSatish Balay   PetscScalar    **rbuf4,*rbuf4_i=NULL,**sbuf_aa,*vals,*mat_a=NULL,*imat_a=NULL,*sbuf_aa_i;
63817df9f7cSHong Zhang   PetscMPIInt    *onodes1,*olengths1,end;
63980d31651SHong Zhang   PetscInt       **row2proc,*row2proc_i,*imat_ilen,*imat_j,*imat_i;
6405c39f6d9SHong Zhang   Mat_SubSppt    *smat_i;
64117df9f7cSHong Zhang   PetscBool      *issorted,colflag,iscsorted=PETSC_TRUE;
64217df9f7cSHong Zhang   PetscInt       *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen;
643a42ab0d6SHong Zhang   PetscInt       bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs;
644a42ab0d6SHong Zhang   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
645a42ab0d6SHong Zhang   PetscInt       nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB;
646afb3cbd8SHong Zhang   PetscScalar    *vworkA=NULL,*vworkB=NULL,*a_a = a->a,*b_a = b->a;
647a42ab0d6SHong Zhang   PetscInt       cstart = c->cstartbs,*bmap = c->garray;
64880d31651SHong Zhang   PetscBool      *allrows,*allcolumns;
649a42ab0d6SHong Zhang 
65017df9f7cSHong Zhang   PetscFunctionBegin;
65117df9f7cSHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
65217df9f7cSHong Zhang   size = c->size;
65317df9f7cSHong Zhang   rank = c->rank;
65417df9f7cSHong Zhang 
6554698041cSHong Zhang   ierr = PetscMalloc5(ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax+1,&allcolumns,ismax,&allrows);CHKERRQ(ierr);
656f489ac74SBarry Smith   ierr = PetscMalloc5(ismax,(PetscInt***)&irow,ismax,(PetscInt***)&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr);
65717df9f7cSHong Zhang 
65817df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
65917df9f7cSHong Zhang     ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr);
66080d31651SHong Zhang     if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */
66117df9f7cSHong Zhang     ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr);
66217df9f7cSHong Zhang 
663ec3c739cSHong Zhang     /* Check for special case: allcolumns */
66417df9f7cSHong Zhang     ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr);
66517df9f7cSHong Zhang     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
6665dd0c08cSHong Zhang 
6671abc2f7bSHong Zhang     if (colflag && ncol[i] == c->Nbs) {
6681abc2f7bSHong Zhang       allcolumns[i] = PETSC_TRUE;
6691abc2f7bSHong Zhang       icol[i]       = NULL;
6701abc2f7bSHong Zhang     } else {
67117df9f7cSHong Zhang       allcolumns[i] = PETSC_FALSE;
67217df9f7cSHong Zhang       ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
6731abc2f7bSHong Zhang     }
674ec3c739cSHong Zhang 
675ec3c739cSHong Zhang     /* Check for special case: allrows */
676ec3c739cSHong Zhang     ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr);
677ec3c739cSHong Zhang     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
678ec3c739cSHong Zhang     if (colflag && nrow[i] == c->Mbs) {
679ec3c739cSHong Zhang       allrows[i] = PETSC_TRUE;
680ec3c739cSHong Zhang       irow[i]    = NULL;
681ec3c739cSHong Zhang     } else {
682ec3c739cSHong Zhang       allrows[i] = PETSC_FALSE;
683ec3c739cSHong Zhang       ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
684ec3c739cSHong Zhang     }
68517df9f7cSHong Zhang   }
68617df9f7cSHong Zhang 
68717df9f7cSHong Zhang   if (scall == MAT_REUSE_MATRIX) {
68817df9f7cSHong Zhang     /* Assumes new rows are same length as the old rows */
68917df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
69017df9f7cSHong Zhang       subc = (Mat_SeqBAIJ*)(submats[i]->data);
6919e686edcSHong Zhang       if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
69217df9f7cSHong Zhang 
69317df9f7cSHong Zhang       /* Initial matrix as if empty */
694580bdb30SBarry Smith       ierr = PetscArrayzero(subc->ilen,subc->mbs);CHKERRQ(ierr);
69517df9f7cSHong Zhang 
69617df9f7cSHong Zhang       /* Initial matrix as if empty */
69717df9f7cSHong Zhang       submats[i]->factortype = C->factortype;
69817df9f7cSHong Zhang 
69917df9f7cSHong Zhang       smat_i   = subc->submatis1;
70017df9f7cSHong Zhang 
70117df9f7cSHong Zhang       nrqs        = smat_i->nrqs;
70217df9f7cSHong Zhang       nrqr        = smat_i->nrqr;
70317df9f7cSHong Zhang       rbuf1       = smat_i->rbuf1;
70417df9f7cSHong Zhang       rbuf2       = smat_i->rbuf2;
70517df9f7cSHong Zhang       rbuf3       = smat_i->rbuf3;
70617df9f7cSHong Zhang       req_source2 = smat_i->req_source2;
70717df9f7cSHong Zhang 
70817df9f7cSHong Zhang       sbuf1     = smat_i->sbuf1;
70917df9f7cSHong Zhang       sbuf2     = smat_i->sbuf2;
71017df9f7cSHong Zhang       ptr       = smat_i->ptr;
71117df9f7cSHong Zhang       tmp       = smat_i->tmp;
71217df9f7cSHong Zhang       ctr       = smat_i->ctr;
71317df9f7cSHong Zhang 
71417df9f7cSHong Zhang       pa          = smat_i->pa;
71517df9f7cSHong Zhang       req_size    = smat_i->req_size;
71617df9f7cSHong Zhang       req_source1 = smat_i->req_source1;
71717df9f7cSHong Zhang 
71817df9f7cSHong Zhang       allcolumns[i] = smat_i->allcolumns;
719ec3c739cSHong Zhang       allrows[i]    = smat_i->allrows;
72017df9f7cSHong Zhang       row2proc[i]   = smat_i->row2proc;
72117df9f7cSHong Zhang       rmap[i]       = smat_i->rmap;
72217df9f7cSHong Zhang       cmap[i]       = smat_i->cmap;
72317df9f7cSHong Zhang     }
7244698041cSHong Zhang 
7254698041cSHong Zhang     if (!ismax){ /* Get dummy submatrices and retrieve struct submatis1 */
7264698041cSHong Zhang       if (!submats[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"submats are null, cannot reuse");
7275c39f6d9SHong Zhang       smat_i = (Mat_SubSppt*)submats[0]->data;
7284698041cSHong Zhang 
7294698041cSHong Zhang       nrqs        = smat_i->nrqs;
7304698041cSHong Zhang       nrqr        = smat_i->nrqr;
7314698041cSHong Zhang       rbuf1       = smat_i->rbuf1;
7324698041cSHong Zhang       rbuf2       = smat_i->rbuf2;
7334698041cSHong Zhang       rbuf3       = smat_i->rbuf3;
7344698041cSHong Zhang       req_source2 = smat_i->req_source2;
7354698041cSHong Zhang 
7364698041cSHong Zhang       sbuf1       = smat_i->sbuf1;
7374698041cSHong Zhang       sbuf2       = smat_i->sbuf2;
7384698041cSHong Zhang       ptr         = smat_i->ptr;
7394698041cSHong Zhang       tmp         = smat_i->tmp;
7404698041cSHong Zhang       ctr         = smat_i->ctr;
7414698041cSHong Zhang 
7424698041cSHong Zhang       pa          = smat_i->pa;
7434698041cSHong Zhang       req_size    = smat_i->req_size;
7444698041cSHong Zhang       req_source1 = smat_i->req_source1;
7454698041cSHong Zhang 
746c67fe082SHong Zhang       allcolumns[0] = PETSC_FALSE;
7474698041cSHong Zhang     }
74817df9f7cSHong Zhang   } else { /* scall == MAT_INITIAL_MATRIX */
74917df9f7cSHong Zhang     /* Get some new tags to keep the communication clean */
75017df9f7cSHong Zhang     ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
75117df9f7cSHong Zhang     ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
75217df9f7cSHong Zhang 
75317df9f7cSHong Zhang     /* evaluate communication - mesg to who, length of mesg, and buffer space
75417df9f7cSHong Zhang      required. Based on this, buffers are allocated, and data copied into them*/
75517df9f7cSHong Zhang     ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr);   /* mesg size, initialize work vectors */
75617df9f7cSHong Zhang 
75717df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
75817df9f7cSHong Zhang       jmax   = nrow[i];
75917df9f7cSHong Zhang       irow_i = irow[i];
76017df9f7cSHong Zhang 
76117df9f7cSHong Zhang       ierr   = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr);
76217df9f7cSHong Zhang       row2proc[i] = row2proc_i;
76317df9f7cSHong Zhang 
76417df9f7cSHong Zhang       if (issorted[i]) proc = 0;
76517df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
76617df9f7cSHong Zhang         if (!issorted[i]) proc = 0;
767ec3c739cSHong Zhang         if (allrows[i]) row = j;
768ec3c739cSHong Zhang         else row = irow_i[j];
769ec3c739cSHong Zhang 
770a42ab0d6SHong Zhang         while (row >= c->rangebs[proc+1]) proc++;
77117df9f7cSHong Zhang         w4[proc]++;
77217df9f7cSHong Zhang         row2proc_i[j] = proc; /* map row index to proc */
77317df9f7cSHong Zhang       }
77417df9f7cSHong Zhang       for (j=0; j<size; j++) {
77517df9f7cSHong Zhang         if (w4[j]) { w1[j] += w4[j];  w3[j]++; w4[j] = 0;}
77617df9f7cSHong Zhang       }
77717df9f7cSHong Zhang     }
77817df9f7cSHong Zhang 
77917df9f7cSHong Zhang     nrqs     = 0;              /* no of outgoing messages */
78017df9f7cSHong Zhang     msz      = 0;              /* total mesg length (for all procs) */
78117df9f7cSHong Zhang     w1[rank] = 0;              /* no mesg sent to self */
78217df9f7cSHong Zhang     w3[rank] = 0;
78317df9f7cSHong Zhang     for (i=0; i<size; i++) {
78417df9f7cSHong Zhang       if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
78517df9f7cSHong Zhang     }
78617df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/
78717df9f7cSHong Zhang     for (i=0,j=0; i<size; i++) {
78817df9f7cSHong Zhang       if (w1[i]) { pa[j] = i; j++; }
78917df9f7cSHong Zhang     }
79017df9f7cSHong Zhang 
79117df9f7cSHong Zhang     /* Each message would have a header = 1 + 2*(no of IS) + data */
79217df9f7cSHong Zhang     for (i=0; i<nrqs; i++) {
79317df9f7cSHong Zhang       j      = pa[i];
79417df9f7cSHong Zhang       w1[j] += w2[j] + 2* w3[j];
79517df9f7cSHong Zhang       msz   += w1[j];
79617df9f7cSHong Zhang     }
79717df9f7cSHong Zhang     ierr = PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);CHKERRQ(ierr);
79817df9f7cSHong Zhang 
79917df9f7cSHong Zhang     /* Determine the number of messages to expect, their lengths, from from-ids */
80017df9f7cSHong Zhang     ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
80117df9f7cSHong Zhang     ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
80217df9f7cSHong Zhang 
80317df9f7cSHong Zhang     /* Now post the Irecvs corresponding to these messages */
80417df9f7cSHong Zhang     tag0 = ((PetscObject)C)->tag;
80517df9f7cSHong Zhang     ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
80617df9f7cSHong Zhang 
80717df9f7cSHong Zhang     ierr = PetscFree(onodes1);CHKERRQ(ierr);
80817df9f7cSHong Zhang     ierr = PetscFree(olengths1);CHKERRQ(ierr);
80917df9f7cSHong Zhang 
81017df9f7cSHong Zhang     /* Allocate Memory for outgoing messages */
81117df9f7cSHong Zhang     ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr);
812580bdb30SBarry Smith     ierr = PetscArrayzero(sbuf1,size);CHKERRQ(ierr);
813580bdb30SBarry Smith     ierr = PetscArrayzero(ptr,size);CHKERRQ(ierr);
81417df9f7cSHong Zhang 
81517df9f7cSHong Zhang     {
81617df9f7cSHong Zhang       PetscInt *iptr = tmp;
81717df9f7cSHong Zhang       k    = 0;
81817df9f7cSHong Zhang       for (i=0; i<nrqs; i++) {
81917df9f7cSHong Zhang         j        = pa[i];
82017df9f7cSHong Zhang         iptr    += k;
82117df9f7cSHong Zhang         sbuf1[j] = iptr;
82217df9f7cSHong Zhang         k        = w1[j];
82317df9f7cSHong Zhang       }
82417df9f7cSHong Zhang     }
82517df9f7cSHong Zhang 
82617df9f7cSHong Zhang     /* Form the outgoing messages. Initialize the header space */
82717df9f7cSHong Zhang     for (i=0; i<nrqs; i++) {
82817df9f7cSHong Zhang       j           = pa[i];
82917df9f7cSHong Zhang       sbuf1[j][0] = 0;
830580bdb30SBarry Smith       ierr        = PetscArrayzero(sbuf1[j]+1,2*w3[j]);CHKERRQ(ierr);
83117df9f7cSHong Zhang       ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
83217df9f7cSHong Zhang     }
83317df9f7cSHong Zhang 
83417df9f7cSHong Zhang     /* Parse the isrow and copy data into outbuf */
83517df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
83617df9f7cSHong Zhang       row2proc_i = row2proc[i];
837580bdb30SBarry Smith       ierr   = PetscArrayzero(ctr,size);CHKERRQ(ierr);
83817df9f7cSHong Zhang       irow_i = irow[i];
83917df9f7cSHong Zhang       jmax   = nrow[i];
84017df9f7cSHong Zhang       for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
84117df9f7cSHong Zhang         proc = row2proc_i[j];
842ec3c739cSHong Zhang         if (allrows[i]) row = j;
843ec3c739cSHong Zhang         else row = irow_i[j];
844ec3c739cSHong Zhang 
84517df9f7cSHong Zhang         if (proc != rank) { /* copy to the outgoing buf*/
84617df9f7cSHong Zhang           ctr[proc]++;
847ec3c739cSHong Zhang           *ptr[proc] = row;
84817df9f7cSHong Zhang           ptr[proc]++;
84917df9f7cSHong Zhang         }
85017df9f7cSHong Zhang       }
85117df9f7cSHong Zhang       /* Update the headers for the current IS */
85217df9f7cSHong Zhang       for (j=0; j<size; j++) { /* Can Optimise this loop too */
85317df9f7cSHong Zhang         if ((ctr_j = ctr[j])) {
85417df9f7cSHong Zhang           sbuf1_j        = sbuf1[j];
85517df9f7cSHong Zhang           k              = ++sbuf1_j[0];
85617df9f7cSHong Zhang           sbuf1_j[2*k]   = ctr_j;
85717df9f7cSHong Zhang           sbuf1_j[2*k-1] = i;
85817df9f7cSHong Zhang         }
85917df9f7cSHong Zhang       }
86017df9f7cSHong Zhang     }
86117df9f7cSHong Zhang 
86217df9f7cSHong Zhang     /*  Now  post the sends */
86317df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr);
86417df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
86517df9f7cSHong Zhang       j    = pa[i];
866ffc4695bSBarry Smith       ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRMPI(ierr);
86717df9f7cSHong Zhang     }
86817df9f7cSHong Zhang 
86917df9f7cSHong Zhang     /* Post Receives to capture the buffer size */
87017df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr);
87117df9f7cSHong Zhang     ierr = PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);CHKERRQ(ierr);
87217df9f7cSHong Zhang     rbuf2[0] = tmp + msz;
87317df9f7cSHong Zhang     for (i=1; i<nrqs; ++i) {
87417df9f7cSHong Zhang       rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
87517df9f7cSHong Zhang     }
87617df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
87717df9f7cSHong Zhang       j    = pa[i];
878ffc4695bSBarry Smith       ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRMPI(ierr);
87917df9f7cSHong Zhang     }
88017df9f7cSHong Zhang 
88117df9f7cSHong Zhang     /* Send to other procs the buf size they should allocate */
88217df9f7cSHong Zhang     /* Receive messages*/
88317df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr);
88417df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr);
88517df9f7cSHong Zhang     ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr);
88617df9f7cSHong Zhang 
887ffc4695bSBarry Smith     ierr = MPI_Waitall(nrqr,r_waits1,r_status1);CHKERRMPI(ierr);
88817df9f7cSHong Zhang     for (i=0; i<nrqr; ++i) {
88917df9f7cSHong Zhang       req_size[i] = 0;
89017df9f7cSHong Zhang       rbuf1_i        = rbuf1[i];
89117df9f7cSHong Zhang       start          = 2*rbuf1_i[0] + 1;
89255b25c41SPierre Jolivet       ierr           = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRMPI(ierr);
89317df9f7cSHong Zhang       ierr           = PetscMalloc1(end+1,&sbuf2[i]);CHKERRQ(ierr);
89417df9f7cSHong Zhang       sbuf2_i        = sbuf2[i];
89517df9f7cSHong Zhang       for (j=start; j<end; j++) {
896ddb3d6bcSHong Zhang         row             = rbuf1_i[j] - rstart;
897ddb3d6bcSHong Zhang         ncols           = a_i[row+1] - a_i[row] + b_i[row+1] - b_i[row];
89817df9f7cSHong Zhang         sbuf2_i[j]      = ncols;
89917df9f7cSHong Zhang         req_size[i] += ncols;
90017df9f7cSHong Zhang       }
90117df9f7cSHong Zhang       req_source1[i] = r_status1[i].MPI_SOURCE;
90217df9f7cSHong Zhang       /* form the header */
90317df9f7cSHong Zhang       sbuf2_i[0] = req_size[i];
90417df9f7cSHong Zhang       for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
90517df9f7cSHong Zhang 
906ffc4695bSBarry Smith       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRMPI(ierr);
90717df9f7cSHong Zhang     }
908ddb3d6bcSHong Zhang 
90917df9f7cSHong Zhang     ierr = PetscFree(r_status1);CHKERRQ(ierr);
91017df9f7cSHong Zhang     ierr = PetscFree(r_waits1);CHKERRQ(ierr);
91117df9f7cSHong Zhang     ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
91217df9f7cSHong Zhang 
91317df9f7cSHong Zhang     /* Receive messages*/
91417df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr);
91517df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr);
91617df9f7cSHong Zhang 
917ffc4695bSBarry Smith     ierr = MPI_Waitall(nrqs,r_waits2,r_status2);CHKERRMPI(ierr);
91817df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
91917df9f7cSHong Zhang       ierr = PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);CHKERRQ(ierr);
92017df9f7cSHong Zhang       req_source2[i] = r_status2[i].MPI_SOURCE;
921ffc4695bSBarry Smith       ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRMPI(ierr);
92217df9f7cSHong Zhang     }
92317df9f7cSHong Zhang     ierr = PetscFree(r_status2);CHKERRQ(ierr);
92417df9f7cSHong Zhang     ierr = PetscFree(r_waits2);CHKERRQ(ierr);
92517df9f7cSHong Zhang 
92617df9f7cSHong Zhang     /* Wait on sends1 and sends2 */
92717df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr);
92817df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr);
92917df9f7cSHong Zhang 
930ffc4695bSBarry Smith     if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRMPI(ierr);}
931ffc4695bSBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRMPI(ierr);}
93217df9f7cSHong Zhang     ierr = PetscFree(s_status1);CHKERRQ(ierr);
93317df9f7cSHong Zhang     ierr = PetscFree(s_status2);CHKERRQ(ierr);
93417df9f7cSHong Zhang     ierr = PetscFree(s_waits1);CHKERRQ(ierr);
93517df9f7cSHong Zhang     ierr = PetscFree(s_waits2);CHKERRQ(ierr);
93617df9f7cSHong Zhang 
93717df9f7cSHong Zhang     /* Now allocate sending buffers for a->j, and send them off */
93817df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr);
93917df9f7cSHong Zhang     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
94017df9f7cSHong Zhang     ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr);
94117df9f7cSHong Zhang     for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
94217df9f7cSHong Zhang 
94317df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr);
94417df9f7cSHong Zhang     {
94517df9f7cSHong Zhang 
94617df9f7cSHong Zhang       for (i=0; i<nrqr; i++) {
94717df9f7cSHong Zhang         rbuf1_i   = rbuf1[i];
94817df9f7cSHong Zhang         sbuf_aj_i = sbuf_aj[i];
94917df9f7cSHong Zhang         ct1       = 2*rbuf1_i[0] + 1;
95017df9f7cSHong Zhang         ct2       = 0;
95117df9f7cSHong Zhang         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
95217df9f7cSHong Zhang           kmax = rbuf1[i][2*j];
95317df9f7cSHong Zhang           for (k=0; k<kmax; k++,ct1++) {
95417df9f7cSHong Zhang             row    = rbuf1_i[ct1] - rstart;
95517df9f7cSHong Zhang             nzA    = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row];
95617df9f7cSHong Zhang             ncols  = nzA + nzB;
95717df9f7cSHong Zhang             cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
95817df9f7cSHong Zhang 
95917df9f7cSHong Zhang             /* load the column indices for this row into cols */
96017df9f7cSHong Zhang             cols = sbuf_aj_i + ct2;
96117df9f7cSHong Zhang             for (l=0; l<nzB; l++) {
962a42ab0d6SHong Zhang               if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
963a42ab0d6SHong Zhang               else break;
96417df9f7cSHong Zhang             }
965a42ab0d6SHong Zhang             imark = l;
966a42ab0d6SHong Zhang             for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];}
967a42ab0d6SHong Zhang             for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
96817df9f7cSHong Zhang             ct2 += ncols;
96917df9f7cSHong Zhang           }
97017df9f7cSHong Zhang         }
971ffc4695bSBarry Smith         ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRMPI(ierr);
97217df9f7cSHong Zhang       }
97317df9f7cSHong Zhang     }
97417df9f7cSHong Zhang     ierr = PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);CHKERRQ(ierr);
97517df9f7cSHong Zhang 
97617df9f7cSHong Zhang     /* create col map: global col of C -> local col of submatrices */
97717df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
97817df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
97917df9f7cSHong Zhang       if (!allcolumns[i]) {
980a42ab0d6SHong Zhang         ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr);
98117df9f7cSHong Zhang 
98217df9f7cSHong Zhang         jmax   = ncol[i];
98317df9f7cSHong Zhang         icol_i = icol[i];
98417df9f7cSHong Zhang         cmap_i = cmap[i];
98517df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
98617df9f7cSHong Zhang           ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
98717df9f7cSHong Zhang         }
98817df9f7cSHong Zhang       } else cmap[i] = NULL;
98917df9f7cSHong Zhang     }
99017df9f7cSHong Zhang #else
99117df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
99217df9f7cSHong Zhang       if (!allcolumns[i]) {
993a42ab0d6SHong Zhang         ierr   = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr);
99417df9f7cSHong Zhang         jmax   = ncol[i];
99517df9f7cSHong Zhang         icol_i = icol[i];
99617df9f7cSHong Zhang         cmap_i = cmap[i];
997a42ab0d6SHong Zhang         for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1;
99817df9f7cSHong Zhang       } else cmap[i] = NULL;
99917df9f7cSHong Zhang     }
100017df9f7cSHong Zhang #endif
100117df9f7cSHong Zhang 
100217df9f7cSHong Zhang     /* Create lens which is required for MatCreate... */
100317df9f7cSHong Zhang     for (i=0,j=0; i<ismax; i++) j += nrow[i];
100417df9f7cSHong Zhang     ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr);
100517df9f7cSHong Zhang 
100617df9f7cSHong Zhang     if (ismax) {
100717df9f7cSHong Zhang       ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr);
100817df9f7cSHong Zhang     }
100917df9f7cSHong Zhang     for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];
101017df9f7cSHong Zhang 
101117df9f7cSHong Zhang     /* Update lens from local data */
101217df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
101317df9f7cSHong Zhang       row2proc_i = row2proc[i];
101417df9f7cSHong Zhang       jmax = nrow[i];
101517df9f7cSHong Zhang       if (!allcolumns[i]) cmap_i = cmap[i];
101617df9f7cSHong Zhang       irow_i = irow[i];
101717df9f7cSHong Zhang       lens_i = lens[i];
101817df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
1019ec3c739cSHong Zhang         if (allrows[i]) row = j;
1020ec3c739cSHong Zhang         else row = irow_i[j]; /* global blocked row of C */
1021ec3c739cSHong Zhang 
102217df9f7cSHong Zhang         proc = row2proc_i[j];
102317df9f7cSHong Zhang         if (proc == rank) {
1024a42ab0d6SHong Zhang           /* Get indices from matA and then from matB */
102517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
1026a42ab0d6SHong Zhang           PetscInt   tt;
1027a42ab0d6SHong Zhang #endif
1028a42ab0d6SHong Zhang           row    = row - rstart;
1029a42ab0d6SHong Zhang           nzA    = a_i[row+1] - a_i[row];
1030a42ab0d6SHong Zhang           nzB    = b_i[row+1] - b_i[row];
1031a42ab0d6SHong Zhang           cworkA =  a_j + a_i[row];
1032a42ab0d6SHong Zhang           cworkB = b_j + b_i[row];
1033a42ab0d6SHong Zhang 
1034a42ab0d6SHong Zhang           if (!allcolumns[i]) {
1035a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1036a42ab0d6SHong Zhang             for (k=0; k<nzA; k++) {
1037a42ab0d6SHong Zhang               ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1038a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1039a42ab0d6SHong Zhang             }
1040a42ab0d6SHong Zhang             for (k=0; k<nzB; k++) {
1041a42ab0d6SHong Zhang               ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1042a42ab0d6SHong Zhang               if (tt) lens_i[j]++;
1043a42ab0d6SHong Zhang             }
1044a42ab0d6SHong Zhang 
1045a42ab0d6SHong Zhang #else
1046a42ab0d6SHong Zhang             for (k=0; k<nzA; k++) {
1047a42ab0d6SHong Zhang               if (cmap_i[cstart + cworkA[k]]) lens_i[j]++;
1048a42ab0d6SHong Zhang             }
1049a42ab0d6SHong Zhang             for (k=0; k<nzB; k++) {
1050a42ab0d6SHong Zhang               if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++;
1051a42ab0d6SHong Zhang             }
1052a42ab0d6SHong Zhang #endif
1053a42ab0d6SHong Zhang           } else { /* allcolumns */
1054a42ab0d6SHong Zhang             lens_i[j] = nzA + nzB;
1055a42ab0d6SHong Zhang           }
105617df9f7cSHong Zhang         }
105717df9f7cSHong Zhang       }
105817df9f7cSHong Zhang     }
105917df9f7cSHong Zhang 
106017df9f7cSHong Zhang     /* Create row map: global row of C -> local row of submatrices */
106117df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
1062059f6b74SHong Zhang       if (!allrows[i]) {
1063059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE)
1064a42ab0d6SHong Zhang         ierr   = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr);
106517df9f7cSHong Zhang         irow_i = irow[i];
106617df9f7cSHong Zhang         jmax   = nrow[i];
106717df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
1068ec3c739cSHong Zhang           if (allrows[i]) {
1069ec3c739cSHong Zhang             ierr = PetscTableAdd(rmap[i],j+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1070ec3c739cSHong Zhang           } else {
107117df9f7cSHong Zhang             ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
107217df9f7cSHong Zhang           }
107317df9f7cSHong Zhang         }
107417df9f7cSHong Zhang #else
1075a42ab0d6SHong Zhang         ierr   = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr);
107617df9f7cSHong Zhang         rmap_i = rmap[i];
107717df9f7cSHong Zhang         irow_i = irow[i];
107817df9f7cSHong Zhang         jmax   = nrow[i];
107917df9f7cSHong Zhang         for (j=0; j<jmax; j++) {
1080ec3c739cSHong Zhang           if (allrows[i]) rmap_i[j] = j;
1081ec3c739cSHong Zhang           else rmap_i[irow_i[j]] = j;
108217df9f7cSHong Zhang         }
108317df9f7cSHong Zhang #endif
1084059f6b74SHong Zhang       } else rmap[i] = NULL;
1085059f6b74SHong Zhang     }
108617df9f7cSHong Zhang 
108717df9f7cSHong Zhang     /* Update lens from offproc data */
108817df9f7cSHong Zhang     {
108917df9f7cSHong Zhang       PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i;
109017df9f7cSHong Zhang 
109155b25c41SPierre Jolivet       ierr    = MPI_Waitall(nrqs,r_waits3,r_status3);CHKERRMPI(ierr);
109217df9f7cSHong Zhang       for (tmp2=0; tmp2<nrqs; tmp2++) {
109317df9f7cSHong Zhang         sbuf1_i = sbuf1[pa[tmp2]];
109417df9f7cSHong Zhang         jmax    = sbuf1_i[0];
109517df9f7cSHong Zhang         ct1     = 2*jmax+1;
109617df9f7cSHong Zhang         ct2     = 0;
109717df9f7cSHong Zhang         rbuf2_i = rbuf2[tmp2];
109817df9f7cSHong Zhang         rbuf3_i = rbuf3[tmp2];
109917df9f7cSHong Zhang         for (j=1; j<=jmax; j++) {
110017df9f7cSHong Zhang           is_no  = sbuf1_i[2*j-1];
110117df9f7cSHong Zhang           max1   = sbuf1_i[2*j];
110217df9f7cSHong Zhang           lens_i = lens[is_no];
110317df9f7cSHong Zhang           if (!allcolumns[is_no]) cmap_i = cmap[is_no];
110417df9f7cSHong Zhang           rmap_i = rmap[is_no];
110517df9f7cSHong Zhang           for (k=0; k<max1; k++,ct1++) {
1106059f6b74SHong Zhang             if (allrows[is_no]) {
1107059f6b74SHong Zhang               row = sbuf1_i[ct1];
1108059f6b74SHong Zhang             } else {
110917df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
111017df9f7cSHong Zhang               ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
111117df9f7cSHong Zhang               row--;
111217df9f7cSHong Zhang               if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
111317df9f7cSHong Zhang #else
111417df9f7cSHong Zhang               row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
111517df9f7cSHong Zhang #endif
1116059f6b74SHong Zhang             }
111717df9f7cSHong Zhang             max2 = rbuf2_i[ct1];
111817df9f7cSHong Zhang             for (l=0; l<max2; l++,ct2++) {
111917df9f7cSHong Zhang               if (!allcolumns[is_no]) {
112017df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
112117df9f7cSHong Zhang                 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
112217df9f7cSHong Zhang #else
112317df9f7cSHong Zhang                 tcol = cmap_i[rbuf3_i[ct2]];
112417df9f7cSHong Zhang #endif
112517df9f7cSHong Zhang                 if (tcol) lens_i[row]++;
112617df9f7cSHong Zhang               } else { /* allcolumns */
112717df9f7cSHong Zhang                 lens_i[row]++; /* lens_i[row] += max2 ? */
112817df9f7cSHong Zhang               }
112917df9f7cSHong Zhang             }
113017df9f7cSHong Zhang           }
113117df9f7cSHong Zhang         }
113217df9f7cSHong Zhang       }
113317df9f7cSHong Zhang     }
113417df9f7cSHong Zhang     ierr = PetscFree(r_waits3);CHKERRQ(ierr);
1135ffc4695bSBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRMPI(ierr);}
113617df9f7cSHong Zhang     ierr = PetscFree2(r_status3,s_status3);CHKERRQ(ierr);
113717df9f7cSHong Zhang     ierr = PetscFree(s_waits3);CHKERRQ(ierr);
113817df9f7cSHong Zhang 
113917df9f7cSHong Zhang     /* Create the submatrices */
114017df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
1141a42ab0d6SHong Zhang       PetscInt bs_tmp;
1142a42ab0d6SHong Zhang       if (ijonly) bs_tmp = 1;
1143a42ab0d6SHong Zhang       else        bs_tmp = bs;
114417df9f7cSHong Zhang 
114517df9f7cSHong Zhang       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
1146a42ab0d6SHong Zhang       ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
114717df9f7cSHong Zhang 
114817df9f7cSHong Zhang       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
1149a42ab0d6SHong Zhang       ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr);
1150a42ab0d6SHong Zhang       ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */
115117df9f7cSHong Zhang 
11525c39f6d9SHong Zhang       /* create struct Mat_SubSppt and attached it to submat */
115317df9f7cSHong Zhang       ierr = PetscNew(&smat_i);CHKERRQ(ierr);
115417df9f7cSHong Zhang       subc = (Mat_SeqBAIJ*)submats[i]->data;
115517df9f7cSHong Zhang       subc->submatis1 = smat_i;
115617df9f7cSHong Zhang 
115717df9f7cSHong Zhang       smat_i->destroy          = submats[i]->ops->destroy;
1158f68bb481SHong Zhang       submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ;
115917df9f7cSHong Zhang       submats[i]->factortype   = C->factortype;
116017df9f7cSHong Zhang 
116117df9f7cSHong Zhang       smat_i->id          = i;
116217df9f7cSHong Zhang       smat_i->nrqs        = nrqs;
116317df9f7cSHong Zhang       smat_i->nrqr        = nrqr;
116417df9f7cSHong Zhang       smat_i->rbuf1       = rbuf1;
116517df9f7cSHong Zhang       smat_i->rbuf2       = rbuf2;
116617df9f7cSHong Zhang       smat_i->rbuf3       = rbuf3;
116717df9f7cSHong Zhang       smat_i->sbuf2       = sbuf2;
116817df9f7cSHong Zhang       smat_i->req_source2 = req_source2;
116917df9f7cSHong Zhang 
117017df9f7cSHong Zhang       smat_i->sbuf1       = sbuf1;
117117df9f7cSHong Zhang       smat_i->ptr         = ptr;
117217df9f7cSHong Zhang       smat_i->tmp         = tmp;
117317df9f7cSHong Zhang       smat_i->ctr         = ctr;
117417df9f7cSHong Zhang 
117517df9f7cSHong Zhang       smat_i->pa           = pa;
117617df9f7cSHong Zhang       smat_i->req_size     = req_size;
117717df9f7cSHong Zhang       smat_i->req_source1  = req_source1;
117817df9f7cSHong Zhang 
117917df9f7cSHong Zhang       smat_i->allcolumns  = allcolumns[i];
1180ec3c739cSHong Zhang       smat_i->allrows     = allrows[i];
118117df9f7cSHong Zhang       smat_i->singleis    = PETSC_FALSE;
118217df9f7cSHong Zhang       smat_i->row2proc    = row2proc[i];
118317df9f7cSHong Zhang       smat_i->rmap        = rmap[i];
118417df9f7cSHong Zhang       smat_i->cmap        = cmap[i];
118517df9f7cSHong Zhang     }
118617df9f7cSHong Zhang 
11874698041cSHong Zhang     if (!ismax) { /* Create dummy submats[0] for reuse struct subc */
11884698041cSHong Zhang       ierr = MatCreate(PETSC_COMM_SELF,&submats[0]);CHKERRQ(ierr);
11894698041cSHong Zhang       ierr = MatSetSizes(submats[0],0,0,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
11904698041cSHong Zhang       ierr = MatSetType(submats[0],MATDUMMY);CHKERRQ(ierr);
11914698041cSHong Zhang 
11925c39f6d9SHong Zhang       /* create struct Mat_SubSppt and attached it to submat */
11934698041cSHong Zhang       ierr = PetscNewLog(submats[0],&smat_i);CHKERRQ(ierr);
11944698041cSHong Zhang       submats[0]->data = (void*)smat_i;
11954698041cSHong Zhang 
11964698041cSHong Zhang       smat_i->destroy          = submats[0]->ops->destroy;
1197f68bb481SHong Zhang       submats[0]->ops->destroy = MatDestroySubMatrix_Dummy;
11984698041cSHong Zhang       submats[0]->factortype   = C->factortype;
11994698041cSHong Zhang 
1200006cef31SHong Zhang       smat_i->id          = 0;
12014698041cSHong Zhang       smat_i->nrqs        = nrqs;
12024698041cSHong Zhang       smat_i->nrqr        = nrqr;
12034698041cSHong Zhang       smat_i->rbuf1       = rbuf1;
12044698041cSHong Zhang       smat_i->rbuf2       = rbuf2;
12054698041cSHong Zhang       smat_i->rbuf3       = rbuf3;
12064698041cSHong Zhang       smat_i->sbuf2       = sbuf2;
12074698041cSHong Zhang       smat_i->req_source2 = req_source2;
12084698041cSHong Zhang 
12094698041cSHong Zhang       smat_i->sbuf1       = sbuf1;
12104698041cSHong Zhang       smat_i->ptr         = ptr;
12114698041cSHong Zhang       smat_i->tmp         = tmp;
12124698041cSHong Zhang       smat_i->ctr         = ctr;
12134698041cSHong Zhang 
12144698041cSHong Zhang       smat_i->pa           = pa;
12154698041cSHong Zhang       smat_i->req_size     = req_size;
12164698041cSHong Zhang       smat_i->req_source1  = req_source1;
12174698041cSHong Zhang 
1218c67fe082SHong Zhang       smat_i->allcolumns  = PETSC_FALSE;
12194698041cSHong Zhang       smat_i->singleis    = PETSC_FALSE;
12204698041cSHong Zhang       smat_i->row2proc    = NULL;
12214698041cSHong Zhang       smat_i->rmap        = NULL;
12224698041cSHong Zhang       smat_i->cmap        = NULL;
12234698041cSHong Zhang     }
12244698041cSHong Zhang 
12254698041cSHong Zhang 
122617df9f7cSHong Zhang     if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);}
122717df9f7cSHong Zhang     ierr = PetscFree(lens);CHKERRQ(ierr);
122817df9f7cSHong Zhang     ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
122917df9f7cSHong Zhang     ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
123017df9f7cSHong Zhang 
123117df9f7cSHong Zhang   } /* endof scall == MAT_INITIAL_MATRIX */
123217df9f7cSHong Zhang 
123317df9f7cSHong Zhang   /* Post recv matrix values */
1234bbc89d27SHong Zhang   if (!ijonly) {
123517df9f7cSHong Zhang     ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr);
123617df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr);
123717df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr);
123817df9f7cSHong Zhang     ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr);
123917df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr);
124017df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
1241a42ab0d6SHong Zhang       ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr);
1242ffc4695bSBarry Smith       ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRMPI(ierr);
124317df9f7cSHong Zhang     }
124417df9f7cSHong Zhang 
124517df9f7cSHong Zhang     /* Allocate sending buffers for a->a, and send them off */
124617df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr);
124717df9f7cSHong Zhang     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
12489e686edcSHong Zhang 
1249a42ab0d6SHong Zhang     ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr);
1250a42ab0d6SHong Zhang     for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
125117df9f7cSHong Zhang 
125217df9f7cSHong Zhang     ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr);
125317df9f7cSHong Zhang 
125417df9f7cSHong Zhang     for (i=0; i<nrqr; i++) {
125517df9f7cSHong Zhang       rbuf1_i   = rbuf1[i];
125617df9f7cSHong Zhang       sbuf_aa_i = sbuf_aa[i];
125717df9f7cSHong Zhang       ct1       = 2*rbuf1_i[0]+1;
125817df9f7cSHong Zhang       ct2       = 0;
125917df9f7cSHong Zhang       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
126017df9f7cSHong Zhang         kmax = rbuf1_i[2*j];
126117df9f7cSHong Zhang         for (k=0; k<kmax; k++,ct1++) {
126217df9f7cSHong Zhang           row    = rbuf1_i[ct1] - rstart;
1263a42ab0d6SHong Zhang           nzA    = a_i[row+1] - a_i[row];
1264a42ab0d6SHong Zhang           nzB    = b_i[row+1] - b_i[row];
126517df9f7cSHong Zhang           ncols  = nzA + nzB;
126617df9f7cSHong Zhang           cworkB = b_j + b_i[row];
1267a42ab0d6SHong Zhang           vworkA = a_a + a_i[row]*bs2;
1268a42ab0d6SHong Zhang           vworkB = b_a + b_i[row]*bs2;
126917df9f7cSHong Zhang 
127017df9f7cSHong Zhang           /* load the column values for this row into vals*/
1271a42ab0d6SHong Zhang           vals = sbuf_aa_i+ct2*bs2;
1272a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1273a42ab0d6SHong Zhang             if ((bmap[cworkB[l]]) < cstart) {
1274580bdb30SBarry Smith               ierr = PetscArraycpy(vals+l*bs2,vworkB+l*bs2,bs2);CHKERRQ(ierr);
1275a42ab0d6SHong Zhang             } else break;
1276a42ab0d6SHong Zhang           }
1277a42ab0d6SHong Zhang           imark = l;
1278a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1279580bdb30SBarry Smith             ierr = PetscArraycpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2);CHKERRQ(ierr);
1280a42ab0d6SHong Zhang           }
1281a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1282580bdb30SBarry Smith             ierr = PetscArraycpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2);CHKERRQ(ierr);
1283a42ab0d6SHong Zhang           }
128417df9f7cSHong Zhang 
128517df9f7cSHong Zhang           ct2 += ncols;
128617df9f7cSHong Zhang         }
128717df9f7cSHong Zhang       }
1288ffc4695bSBarry Smith       ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRMPI(ierr);
128917df9f7cSHong Zhang     }
1290bbc89d27SHong Zhang   }
129117df9f7cSHong Zhang 
129217df9f7cSHong Zhang   /* Assemble the matrices */
129317df9f7cSHong Zhang   /* First assemble the local rows */
129417df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
129517df9f7cSHong Zhang     row2proc_i = row2proc[i];
129617df9f7cSHong Zhang     subc      = (Mat_SeqBAIJ*)submats[i]->data;
129717df9f7cSHong Zhang     imat_ilen = subc->ilen;
129817df9f7cSHong Zhang     imat_j    = subc->j;
129917df9f7cSHong Zhang     imat_i    = subc->i;
130017df9f7cSHong Zhang     imat_a    = subc->a;
130117df9f7cSHong Zhang 
130217df9f7cSHong Zhang     if (!allcolumns[i]) cmap_i = cmap[i];
130317df9f7cSHong Zhang     rmap_i = rmap[i];
130417df9f7cSHong Zhang     irow_i = irow[i];
130517df9f7cSHong Zhang     jmax   = nrow[i];
130617df9f7cSHong Zhang     for (j=0; j<jmax; j++) {
1307ec3c739cSHong Zhang       if (allrows[i]) row = j;
1308ec3c739cSHong Zhang       else row  = irow_i[j];
130917df9f7cSHong Zhang       proc = row2proc_i[j];
1310a42ab0d6SHong Zhang 
131117df9f7cSHong Zhang       if (proc == rank) {
1312a42ab0d6SHong Zhang 
1313a42ab0d6SHong Zhang         row    = row - rstart;
1314a42ab0d6SHong Zhang         nzA    = a_i[row+1] - a_i[row];
1315a42ab0d6SHong Zhang         nzB    = b_i[row+1] - b_i[row];
1316a42ab0d6SHong Zhang         cworkA = a_j + a_i[row];
1317a42ab0d6SHong Zhang         cworkB = b_j + b_i[row];
1318a42ab0d6SHong Zhang         if (!ijonly) {
1319a42ab0d6SHong Zhang           vworkA = a_a + a_i[row]*bs2;
1320a42ab0d6SHong Zhang           vworkB = b_a + b_i[row]*bs2;
1321a42ab0d6SHong Zhang         }
1322059f6b74SHong Zhang 
1323059f6b74SHong Zhang         if (allrows[i]) {
1324059f6b74SHong Zhang           row = row+rstart;
1325059f6b74SHong Zhang         } else {
1326a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1327a42ab0d6SHong Zhang           ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr);
1328a42ab0d6SHong Zhang           row--;
1329121971b2SHong Zhang 
1330a42ab0d6SHong Zhang           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1331a42ab0d6SHong Zhang #else
1332a42ab0d6SHong Zhang           row = rmap_i[row + rstart];
1333a42ab0d6SHong Zhang #endif
1334059f6b74SHong Zhang         }
1335a42ab0d6SHong Zhang         mat_i = imat_i[row];
1336a42ab0d6SHong Zhang         if (!ijonly) mat_a = imat_a + mat_i*bs2;
1337a42ab0d6SHong Zhang         mat_j    = imat_j + mat_i;
133880d31651SHong Zhang         ilen = imat_ilen[row];
1339a42ab0d6SHong Zhang 
1340a42ab0d6SHong Zhang         /* load the column indices for this row into cols*/
1341a42ab0d6SHong Zhang         if (!allcolumns[i]) {
1342a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1343a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1344a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1345a42ab0d6SHong Zhang               ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr);
1346a42ab0d6SHong Zhang               if (tcol) {
1347a42ab0d6SHong Zhang #else
1348a42ab0d6SHong Zhang               if ((tcol = cmap_i[ctmp])) {
1349a42ab0d6SHong Zhang #endif
1350a42ab0d6SHong Zhang                 *mat_j++ = tcol - 1;
1351580bdb30SBarry Smith                 ierr     = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr);
1352a42ab0d6SHong Zhang                 mat_a   += bs2;
135380d31651SHong Zhang                 ilen++;
1354a42ab0d6SHong Zhang               }
1355a42ab0d6SHong Zhang             } else break;
1356a42ab0d6SHong Zhang           }
1357a42ab0d6SHong Zhang           imark = l;
1358a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1359a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1360a42ab0d6SHong Zhang             ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1361a42ab0d6SHong Zhang             if (tcol) {
1362a42ab0d6SHong Zhang #else
1363a42ab0d6SHong Zhang             if ((tcol = cmap_i[cstart + cworkA[l]])) {
1364a42ab0d6SHong Zhang #endif
1365a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1366a42ab0d6SHong Zhang               if (!ijonly) {
1367580bdb30SBarry Smith                 ierr   = PetscArraycpy(mat_a,vworkA+l*bs2,bs2);CHKERRQ(ierr);
1368a42ab0d6SHong Zhang                 mat_a += bs2;
1369a42ab0d6SHong Zhang               }
137080d31651SHong Zhang               ilen++;
1371a42ab0d6SHong Zhang             }
1372a42ab0d6SHong Zhang           }
1373a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1374a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE)
1375a42ab0d6SHong Zhang             ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1376a42ab0d6SHong Zhang             if (tcol) {
1377a42ab0d6SHong Zhang #else
1378a42ab0d6SHong Zhang             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1379a42ab0d6SHong Zhang #endif
1380a42ab0d6SHong Zhang               *mat_j++ = tcol - 1;
1381a42ab0d6SHong Zhang               if (!ijonly) {
1382580bdb30SBarry Smith                 ierr   = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr);
1383a42ab0d6SHong Zhang                 mat_a += bs2;
1384a42ab0d6SHong Zhang               }
138580d31651SHong Zhang               ilen++;
1386a42ab0d6SHong Zhang             }
1387a42ab0d6SHong Zhang           }
1388a42ab0d6SHong Zhang         } else { /* allcolumns */
1389a42ab0d6SHong Zhang           for (l=0; l<nzB; l++) {
1390a42ab0d6SHong Zhang             if ((ctmp = bmap[cworkB[l]]) < cstart) {
1391a42ab0d6SHong Zhang               *mat_j++ = ctmp;
1392580bdb30SBarry Smith               ierr     = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr);
1393a42ab0d6SHong Zhang               mat_a   += bs2;
139480d31651SHong Zhang               ilen++;
1395a42ab0d6SHong Zhang             } else break;
1396a42ab0d6SHong Zhang           }
1397a42ab0d6SHong Zhang           imark = l;
1398a42ab0d6SHong Zhang           for (l=0; l<nzA; l++) {
1399a42ab0d6SHong Zhang             *mat_j++ = cstart+cworkA[l];
1400a42ab0d6SHong Zhang             if (!ijonly) {
1401580bdb30SBarry Smith               ierr   = PetscArraycpy(mat_a,vworkA+l*bs2,bs2);CHKERRQ(ierr);
1402a42ab0d6SHong Zhang               mat_a += bs2;
1403a42ab0d6SHong Zhang             }
140480d31651SHong Zhang             ilen++;
1405a42ab0d6SHong Zhang           }
1406a42ab0d6SHong Zhang           for (l=imark; l<nzB; l++) {
1407a42ab0d6SHong Zhang             *mat_j++ = bmap[cworkB[l]];
1408a42ab0d6SHong Zhang             if (!ijonly) {
1409580bdb30SBarry Smith               ierr   = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr);
1410a42ab0d6SHong Zhang               mat_a += bs2;
1411a42ab0d6SHong Zhang             }
141280d31651SHong Zhang             ilen++;
1413a42ab0d6SHong Zhang           }
1414a42ab0d6SHong Zhang         }
141580d31651SHong Zhang         imat_ilen[row] = ilen;
141617df9f7cSHong Zhang       }
141717df9f7cSHong Zhang     }
141817df9f7cSHong Zhang   }
141917df9f7cSHong Zhang 
142017df9f7cSHong Zhang   /* Now assemble the off proc rows */
14215dd0c08cSHong Zhang   if (!ijonly) {
1422ffc4695bSBarry Smith     ierr = MPI_Waitall(nrqs,r_waits4,r_status4);CHKERRMPI(ierr);
14235dd0c08cSHong Zhang   }
142417df9f7cSHong Zhang   for (tmp2=0; tmp2<nrqs; tmp2++) {
142517df9f7cSHong Zhang     sbuf1_i = sbuf1[pa[tmp2]];
142617df9f7cSHong Zhang     jmax    = sbuf1_i[0];
142717df9f7cSHong Zhang     ct1     = 2*jmax + 1;
142817df9f7cSHong Zhang     ct2     = 0;
142917df9f7cSHong Zhang     rbuf2_i = rbuf2[tmp2];
143017df9f7cSHong Zhang     rbuf3_i = rbuf3[tmp2];
14315dd0c08cSHong Zhang     if (!ijonly) rbuf4_i = rbuf4[tmp2];
143217df9f7cSHong Zhang     for (j=1; j<=jmax; j++) {
143317df9f7cSHong Zhang       is_no     = sbuf1_i[2*j-1];
143417df9f7cSHong Zhang       rmap_i    = rmap[is_no];
143517df9f7cSHong Zhang       if (!allcolumns[is_no]) cmap_i = cmap[is_no];
143617df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ*)submats[is_no]->data;
143717df9f7cSHong Zhang       imat_ilen = subc->ilen;
143817df9f7cSHong Zhang       imat_j    = subc->j;
143917df9f7cSHong Zhang       imat_i    = subc->i;
14405dd0c08cSHong Zhang       if (!ijonly) imat_a    = subc->a;
144117df9f7cSHong Zhang       max1      = sbuf1_i[2*j];
14425dd0c08cSHong Zhang       for (k=0; k<max1; k++,ct1++) { /* for each recved block row */
144317df9f7cSHong Zhang         row = sbuf1_i[ct1];
1444059f6b74SHong Zhang 
1445059f6b74SHong Zhang         if (allrows[is_no]) {
1446059f6b74SHong Zhang           row = sbuf1_i[ct1];
1447059f6b74SHong Zhang         } else {
144817df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
144917df9f7cSHong Zhang           ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr);
145017df9f7cSHong Zhang           row--;
14515dd0c08cSHong Zhang           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
145217df9f7cSHong Zhang #else
145317df9f7cSHong Zhang           row = rmap_i[row];
145417df9f7cSHong Zhang #endif
1455059f6b74SHong Zhang         }
145617df9f7cSHong Zhang         ilen  = imat_ilen[row];
145717df9f7cSHong Zhang         mat_i = imat_i[row];
14585dd0c08cSHong Zhang         if (!ijonly) mat_a = imat_a + mat_i*bs2;
145917df9f7cSHong Zhang         mat_j = imat_j + mat_i;
146017df9f7cSHong Zhang         max2  = rbuf2_i[ct1];
146117df9f7cSHong Zhang         if (!allcolumns[is_no]) {
146217df9f7cSHong Zhang           for (l=0; l<max2; l++,ct2++) {
146317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE)
146417df9f7cSHong Zhang             ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
146517df9f7cSHong Zhang #else
146617df9f7cSHong Zhang             tcol = cmap_i[rbuf3_i[ct2]];
146717df9f7cSHong Zhang #endif
146817df9f7cSHong Zhang             if (tcol) {
146917df9f7cSHong Zhang               *mat_j++ = tcol - 1;
14705dd0c08cSHong Zhang               if (!ijonly) {
1471580bdb30SBarry Smith                 ierr   = PetscArraycpy(mat_a,rbuf4_i+ct2*bs2,bs2);CHKERRQ(ierr);
14725dd0c08cSHong Zhang                 mat_a += bs2;
14735dd0c08cSHong Zhang               }
147417df9f7cSHong Zhang               ilen++;
147517df9f7cSHong Zhang             }
147617df9f7cSHong Zhang           }
147717df9f7cSHong Zhang         } else { /* allcolumns */
147817df9f7cSHong Zhang           for (l=0; l<max2; l++,ct2++) {
147917df9f7cSHong Zhang             *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */
14805dd0c08cSHong Zhang             if (!ijonly) {
1481580bdb30SBarry Smith               ierr   = PetscArraycpy(mat_a,rbuf4_i+ct2*bs2,bs2);CHKERRQ(ierr);
14825dd0c08cSHong Zhang               mat_a += bs2;
14835dd0c08cSHong Zhang             }
148417df9f7cSHong Zhang             ilen++;
148517df9f7cSHong Zhang           }
148617df9f7cSHong Zhang         }
148717df9f7cSHong Zhang         imat_ilen[row] = ilen;
148817df9f7cSHong Zhang       }
148917df9f7cSHong Zhang     }
149017df9f7cSHong Zhang   }
149117df9f7cSHong Zhang 
149280d31651SHong Zhang   if (!iscsorted) { /* sort column indices of the rows */
149380d31651SHong Zhang     MatScalar *work;
149480d31651SHong Zhang 
149580d31651SHong Zhang     ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr);
149617df9f7cSHong Zhang     for (i=0; i<ismax; i++) {
149717df9f7cSHong Zhang       subc      = (Mat_SeqBAIJ*)submats[i]->data;
149880d31651SHong Zhang       imat_ilen = subc->ilen;
149917df9f7cSHong Zhang       imat_j    = subc->j;
150017df9f7cSHong Zhang       imat_i    = subc->i;
15014b6ceb0dSHong Zhang       if (!ijonly) imat_a = subc->a;
150217df9f7cSHong Zhang       if (allcolumns[i]) continue;
15034b6ceb0dSHong Zhang 
150417df9f7cSHong Zhang       jmax = nrow[i];
150517df9f7cSHong Zhang       for (j=0; j<jmax; j++) {
150680d31651SHong Zhang         mat_i = imat_i[j];
150780d31651SHong Zhang         mat_j = imat_j + mat_i;
15084b6ceb0dSHong Zhang         ilen  = imat_ilen[j];
150980d31651SHong Zhang         if (ijonly) {
151080d31651SHong Zhang           ierr = PetscSortInt(ilen,mat_j);CHKERRQ(ierr);
151180d31651SHong Zhang         } else {
15124b6ceb0dSHong Zhang           mat_a = imat_a + mat_i*bs2;
151380d31651SHong Zhang           ierr = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);
151417df9f7cSHong Zhang         }
151517df9f7cSHong Zhang       }
151617df9f7cSHong Zhang     }
151780d31651SHong Zhang     ierr = PetscFree(work);CHKERRQ(ierr);
151880d31651SHong Zhang   }
151917df9f7cSHong Zhang 
1520bbc89d27SHong Zhang   if (!ijonly) {
152117df9f7cSHong Zhang     ierr = PetscFree(r_status4);CHKERRQ(ierr);
152217df9f7cSHong Zhang     ierr = PetscFree(r_waits4);CHKERRQ(ierr);
1523ffc4695bSBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRMPI(ierr);}
152417df9f7cSHong Zhang     ierr = PetscFree(s_waits4);CHKERRQ(ierr);
152517df9f7cSHong Zhang     ierr = PetscFree(s_status4);CHKERRQ(ierr);
1526bbc89d27SHong Zhang   }
152717df9f7cSHong Zhang 
152817df9f7cSHong Zhang   /* Restore the indices */
152917df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
1530ec3c739cSHong Zhang     if (!allrows[i]) {
153117df9f7cSHong Zhang       ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1532ec3c739cSHong Zhang     }
153317df9f7cSHong Zhang     if (!allcolumns[i]) {
153417df9f7cSHong Zhang       ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
153517df9f7cSHong Zhang     }
153617df9f7cSHong Zhang   }
153717df9f7cSHong Zhang 
153817df9f7cSHong Zhang   for (i=0; i<ismax; i++) {
153917df9f7cSHong Zhang     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154017df9f7cSHong Zhang     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154117df9f7cSHong Zhang   }
154217df9f7cSHong Zhang 
15432cd485c2SBarry Smith   ierr = PetscFree5(*(PetscInt***)&irow,*(PetscInt***)&icol,nrow,ncol,issorted);CHKERRQ(ierr);
15444698041cSHong Zhang   ierr = PetscFree5(row2proc,cmap,rmap,allcolumns,allrows);CHKERRQ(ierr);
154517df9f7cSHong Zhang 
1546bbc89d27SHong Zhang   if (!ijonly) {
154717df9f7cSHong Zhang     ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
154817df9f7cSHong Zhang     ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
154917df9f7cSHong Zhang 
155017df9f7cSHong Zhang     for (i=0; i<nrqs; ++i) {
155117df9f7cSHong Zhang       ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);
155217df9f7cSHong Zhang     }
155317df9f7cSHong Zhang     ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1554bbc89d27SHong Zhang   }
15557a868f3eSHong Zhang   c->ijonly = PETSC_FALSE; /* set back to the default */
15563a40ed3dSBarry Smith   PetscFunctionReturn(0);
1557d5b485abSSatish Balay }
1558