xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 1795a4d16c893ec2fc06bbbc6c5ce592a2de75d4)
123bdbc58SBarry Smith 
2d5b485abSSatish Balay /*
3d5b485abSSatish Balay    Routines to compute overlapping regions of a parallel MPI matrix
4d5b485abSSatish Balay   and to find submatrices that were shared across processors.
5d5b485abSSatish Balay */
6c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h>
7c6db04a5SJed Brown #include <petscbt.h>
8d5b485abSSatish Balay 
9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char**,PetscInt*,PetscInt**);
10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt**,PetscInt**,PetscInt*);
1109573ac7SBarry Smith extern PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
1209573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
13d5b485abSSatish Balay 
144a2ae208SSatish Balay #undef __FUNCT__
154a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ"
16b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov)
17d5b485abSSatish Balay {
186849ba73SBarry Smith   PetscErrorCode ierr;
19d0f46423SBarry Smith   PetscInt       i,N=C->cmap->N, bs=C->rmap->bs;
2036f4e84dSSatish Balay   IS             *is_new;
2136f4e84dSSatish Balay 
223a40ed3dSBarry Smith   PetscFunctionBegin;
23785e854fSJed Brown   ierr = PetscMalloc1(imax,&is_new);CHKERRQ(ierr);
24df36731bSBarry Smith   /* Convert the indices into block format */
2505d8c843SHong Zhang   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr);
2626fbe8dcSKarl Rupp   if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");
27d5b485abSSatish Balay   for (i=0; i<ov; ++i) {
2836f4e84dSSatish Balay     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr);
29d5b485abSSatish Balay   }
306bf464f9SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);}
3105d8c843SHong Zhang   ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr);
326bf464f9SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);}
33606d414cSSatish Balay   ierr = PetscFree(is_new);CHKERRQ(ierr);
343a40ed3dSBarry Smith   PetscFunctionReturn(0);
35d5b485abSSatish Balay }
36d5b485abSSatish Balay 
37d5b485abSSatish Balay /*
38d5b485abSSatish Balay   Sample message format:
39d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
400e9b0e7eSHong Zhang   to index sets is[1], is[5]
41d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
42d5b485abSSatish Balay   -----------
43d5b485abSSatish Balay   mesg [1] = 1 => is[1]
44d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
45d5b485abSSatish Balay   -----------
46d5b485abSSatish Balay   mesg [5] = 5  => is[5]
47d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
48d5b485abSSatish Balay   -----------
49d5b485abSSatish Balay   mesg [7]
500e9b0e7eSHong Zhang   mesg [n]  data(is[1])
51d5b485abSSatish Balay   -----------
52d5b485abSSatish Balay   mesg[n+1]
53d5b485abSSatish Balay   mesg[m]  data(is[5])
54d5b485abSSatish Balay   -----------
55d5b485abSSatish Balay 
56d5b485abSSatish Balay   Notes:
57d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
58d5b485abSSatish Balay   nrqr - no of requests recieved (which have to be or which have been processed
59d5b485abSSatish Balay */
604a2ae208SSatish Balay #undef __FUNCT__
614a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once"
62db41cccfSHong Zhang PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[])
63d5b485abSSatish Balay {
64df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
655d0c19d7SBarry Smith   const PetscInt **idx,*idx_i;
6624c049a4SHong Zhang   PetscInt       *n,*w3,*w4,**data,len;
676849ba73SBarry Smith   PetscErrorCode ierr;
68b24ad042SBarry Smith   PetscMPIInt    size,rank,tag1,tag2,*w2,*w1,nrqr;
69245d216aSHong Zhang   PetscInt       Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr;
708f9f447aSHong Zhang   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p;
71b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1,*onodes2,*olengths2;
72f1af5d2fSBarry Smith   PetscBT        *table;
73d5b485abSSatish Balay   MPI_Comm       comm;
74d5b485abSSatish Balay   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
75d5b485abSSatish Balay   MPI_Status     *s_status,*recv_status;
768f9f447aSHong Zhang   char           *t_p;
77d5b485abSSatish Balay 
783a40ed3dSBarry Smith   PetscFunctionBegin;
79ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
80d5b485abSSatish Balay   size = c->size;
81d5b485abSSatish Balay   rank = c->rank;
82df36731bSBarry Smith   Mbs  = c->Mbs;
83d5b485abSSatish Balay 
84c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
85c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
86c7dd2924SSatish Balay 
87dcca6d9dSJed Brown   ierr = PetscMalloc2(imax+1,&idx,imax,&n);CHKERRQ(ierr);
8824c049a4SHong Zhang 
89d5b485abSSatish Balay   for (i=0; i<imax; i++) {
90d5b485abSSatish Balay     ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr);
91b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
92d5b485abSSatish Balay   }
93d5b485abSSatish Balay 
94d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
95d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
96dcca6d9dSJed Brown   ierr = PetscMalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr);
9723bdbc58SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
9823bdbc58SBarry Smith   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
9923bdbc58SBarry Smith   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
100d5b485abSSatish Balay   for (i=0; i<imax; i++) {
101b24ad042SBarry Smith     ierr  = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
102d5b485abSSatish Balay     idx_i = idx[i];
103d5b485abSSatish Balay     len   = n[i];
104d5b485abSSatish Balay     for (j=0; j<len; j++) {
105d5b485abSSatish Balay       row = idx_i[j];
106f23aa3ddSBarry Smith       if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries");
107ca31476aSJed Brown       ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr);
108d5b485abSSatish Balay       w4[proc]++;
109d5b485abSSatish Balay     }
110d5b485abSSatish Balay     for (j=0; j<size; j++) {
111d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
112d5b485abSSatish Balay     }
113d5b485abSSatish Balay   }
114d5b485abSSatish Balay 
115d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
116d5b485abSSatish Balay   msz      = 0;              /* total mesg length (for all proc */
1170e9b0e7eSHong Zhang   w1[rank] = 0;              /* no mesg sent to itself */
118d5b485abSSatish Balay   w3[rank] = 0;
119d5b485abSSatish Balay   for (i=0; i<size; i++) {
120d5b485abSSatish Balay     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
121d5b485abSSatish Balay   }
122d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
123785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&pa);CHKERRQ(ierr);
124d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
125d5b485abSSatish Balay     if (w1[i]) {pa[j] = i; j++;}
126d5b485abSSatish Balay   }
127d5b485abSSatish Balay 
128d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
129d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
130d5b485abSSatish Balay     j      = pa[i];
131d5b485abSSatish Balay     w1[j] += w2[j] + 2*w3[j];
132d5b485abSSatish Balay     msz   += w1[j];
133d5b485abSSatish Balay   }
134d5b485abSSatish Balay 
135c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
136a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
137c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
138d5b485abSSatish Balay 
139c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
140c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr);
141d5b485abSSatish Balay 
142d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
143dcca6d9dSJed Brown   ierr = PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr);CHKERRQ(ierr);
14423bdbc58SBarry Smith   ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr);
14523bdbc58SBarry Smith   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
146d5b485abSSatish Balay   {
147b24ad042SBarry Smith     PetscInt *iptr = tmp,ict  = 0;
148d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
149d5b485abSSatish Balay       j         = pa[i];
150d5b485abSSatish Balay       iptr     +=  ict;
151d5b485abSSatish Balay       outdat[j] = iptr;
152d5b485abSSatish Balay       ict       = w1[j];
153d5b485abSSatish Balay     }
154d5b485abSSatish Balay   }
155d5b485abSSatish Balay 
156d5b485abSSatish Balay   /* Form the outgoing messages */
157d5b485abSSatish Balay   /*plug in the headers*/
158d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
159d5b485abSSatish Balay     j            = pa[i];
160d5b485abSSatish Balay     outdat[j][0] = 0;
161b24ad042SBarry Smith     ierr         = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
162d5b485abSSatish Balay     ptr[j]       = outdat[j] + 2*w3[j] + 1;
163d5b485abSSatish Balay   }
164d5b485abSSatish Balay 
165d5b485abSSatish Balay   /* Memory for doing local proc's work*/
166d5b485abSSatish Balay   {
167*1795a4d1SJed Brown     ierr = PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mbs*imax,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,&t_p);CHKERRQ(ierr);
168d5b485abSSatish Balay 
169bbd702dbSSatish Balay     for (i=0; i<imax; i++) {
170b6410449SSatish Balay       table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
171bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
172d5b485abSSatish Balay     }
173d5b485abSSatish Balay   }
174d5b485abSSatish Balay 
175d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
176d5b485abSSatish Balay   {
177b24ad042SBarry Smith     PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j;
178f1af5d2fSBarry Smith     PetscBT  table_i;
179d5b485abSSatish Balay 
180d5b485abSSatish Balay     for (i=0; i<imax; i++) {
181b24ad042SBarry Smith       ierr    = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
182d5b485abSSatish Balay       n_i     = n[i];
183d5b485abSSatish Balay       table_i = table[i];
184d5b485abSSatish Balay       idx_i   = idx[i];
185d5b485abSSatish Balay       data_i  = data[i];
186d5b485abSSatish Balay       isz_i   = isz[i];
187d5b485abSSatish Balay       for (j=0; j<n_i; j++) {   /* parse the indices of each IS */
188d5b485abSSatish Balay         row  = idx_i[j];
189ca31476aSJed Brown         ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr);
190d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
191d5b485abSSatish Balay           ctr[proc]++;
192d5b485abSSatish Balay           *ptr[proc] = row;
193d5b485abSSatish Balay           ptr[proc]++;
194d6b45a43SBarry Smith         } else { /* Update the local table */
19526fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row;
196d5b485abSSatish Balay         }
197d5b485abSSatish Balay       }
198d5b485abSSatish Balay       /* Update the headers for the current IS */
199d5b485abSSatish Balay       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
200d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
201d5b485abSSatish Balay           outdat_j        = outdat[j];
202d5b485abSSatish Balay           k               = ++outdat_j[0];
203d5b485abSSatish Balay           outdat_j[2*k]   = ctr_j;
204d5b485abSSatish Balay           outdat_j[2*k-1] = i;
205d5b485abSSatish Balay         }
206d5b485abSSatish Balay       }
207d5b485abSSatish Balay       isz[i] = isz_i;
208d5b485abSSatish Balay     }
209d5b485abSSatish Balay   }
210d5b485abSSatish Balay 
211d5b485abSSatish Balay   /*  Now  post the sends */
212785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&s_waits1);CHKERRQ(ierr);
213d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
214d5b485abSSatish Balay     j    = pa[i];
215b24ad042SBarry Smith     ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr);
216d5b485abSSatish Balay   }
217d5b485abSSatish Balay 
218d5b485abSSatish Balay   /* No longer need the original indices*/
219d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
220d5b485abSSatish Balay     ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr);
221d5b485abSSatish Balay   }
22224c049a4SHong Zhang   ierr = PetscFree2(idx,n);CHKERRQ(ierr);
223d5b485abSSatish Balay 
224d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
2256bf464f9SBarry Smith     ierr = ISDestroy(&is[i]);CHKERRQ(ierr);
226d5b485abSSatish Balay   }
227d5b485abSSatish Balay 
228d5b485abSSatish Balay   /* Do Local work*/
229df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
230d5b485abSSatish Balay 
231d5b485abSSatish Balay   /* Receive messages*/
232785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&recv_status);CHKERRQ(ierr);
2330c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);}
234d5b485abSSatish Balay 
235785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&s_status);CHKERRQ(ierr);
2360c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);}
237d5b485abSSatish Balay 
238d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
23923bdbc58SBarry Smith   ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr);
24023bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
241d5b485abSSatish Balay 
242785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&xdata);CHKERRQ(ierr);
243785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&isz1);CHKERRQ(ierr);
244df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
245c05d87d6SBarry Smith   ierr = PetscFree(rbuf[0]);CHKERRQ(ierr);
246606d414cSSatish Balay   ierr = PetscFree(rbuf);CHKERRQ(ierr);
247d5b485abSSatish Balay 
248d5b485abSSatish Balay   /* Send the data back*/
249d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
250d5b485abSSatish Balay   {
251b24ad042SBarry Smith     PetscMPIInt *rw1;
252d5b485abSSatish Balay 
253785e854fSJed Brown     ierr = PetscMalloc1(size,&rw1);CHKERRQ(ierr);
254b24ad042SBarry Smith     ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr);
255c7dd2924SSatish Balay 
256d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
257d5b485abSSatish Balay       proc = recv_status[i].MPI_SOURCE;
258e32f2f54SBarry Smith       if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch");
259d5b485abSSatish Balay       rw1[proc] = isz1[i];
260d5b485abSSatish Balay     }
261d5b485abSSatish Balay 
262c7dd2924SSatish Balay     ierr = PetscFree(onodes1);CHKERRQ(ierr);
263c7dd2924SSatish Balay     ierr = PetscFree(olengths1);CHKERRQ(ierr);
264d5b485abSSatish Balay 
265c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
266c7dd2924SSatish Balay     ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr);
26703512dc5SSatish Balay     ierr = PetscFree(rw1);CHKERRQ(ierr);
268c7dd2924SSatish Balay   }
269c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
270c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr);
271d5b485abSSatish Balay 
272d5b485abSSatish Balay   /*  Now  post the sends */
273785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&s_waits2);CHKERRQ(ierr);
274d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
275d5b485abSSatish Balay     j    = recv_status[i].MPI_SOURCE;
276b24ad042SBarry Smith     ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr);
277d5b485abSSatish Balay   }
278d5b485abSSatish Balay 
279d5b485abSSatish Balay   /* receive work done on other processors*/
280d5b485abSSatish Balay   {
281b24ad042SBarry Smith     PetscMPIInt idex;
282b24ad042SBarry Smith     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
283f1af5d2fSBarry Smith     PetscBT     table_i;
284d5b485abSSatish Balay     MPI_Status  *status2;
285d5b485abSSatish Balay 
286169449a0SSatish Balay     ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr);
287d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
288999d9058SBarry Smith       ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr);
289d5b485abSSatish Balay       /* Process the message*/
290999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
291d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
292999d9058SBarry Smith       jmax    = rbuf2[idex][0];
293d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
294d5b485abSSatish Balay         max     = rbuf2_i[2*j];
295d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
296d5b485abSSatish Balay         isz_i   = isz[is_no];
297d5b485abSSatish Balay         data_i  = data[is_no];
298d5b485abSSatish Balay         table_i = table[is_no];
299d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
300d5b485abSSatish Balay           row = rbuf2_i[ct1];
30126fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row;
302d5b485abSSatish Balay         }
303d5b485abSSatish Balay         isz[is_no] = isz_i;
304d5b485abSSatish Balay       }
305d5b485abSSatish Balay     }
3060c468ba9SBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);}
307606d414cSSatish Balay     ierr = PetscFree(status2);CHKERRQ(ierr);
308d5b485abSSatish Balay   }
309d5b485abSSatish Balay 
310d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
31170b3c8c7SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
312d5b485abSSatish Balay   }
313d5b485abSSatish Balay 
314c7dd2924SSatish Balay 
315c7dd2924SSatish Balay   ierr = PetscFree(onodes2);CHKERRQ(ierr);
316c7dd2924SSatish Balay   ierr = PetscFree(olengths2);CHKERRQ(ierr);
317c7dd2924SSatish Balay 
318606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
319c05d87d6SBarry Smith   ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr);
320606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
321606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
322606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
323606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
324606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
3258f9f447aSHong Zhang   ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr);
326606d414cSSatish Balay   ierr = PetscFree(s_status);CHKERRQ(ierr);
327606d414cSSatish Balay   ierr = PetscFree(recv_status);CHKERRQ(ierr);
328606d414cSSatish Balay   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
329606d414cSSatish Balay   ierr = PetscFree(xdata);CHKERRQ(ierr);
330606d414cSSatish Balay   ierr = PetscFree(isz1);CHKERRQ(ierr);
3313a40ed3dSBarry Smith   PetscFunctionReturn(0);
332d5b485abSSatish Balay }
333d5b485abSSatish Balay 
3344a2ae208SSatish Balay #undef __FUNCT__
3354a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local"
336d5b485abSSatish Balay /*
337df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
338d5b485abSSatish Balay        the work on the local processor.
339d5b485abSSatish Balay 
340d5b485abSSatish Balay      Inputs:
341df36731bSBarry Smith       C      - MAT_MPIBAIJ;
342d5b485abSSatish Balay       imax - total no of index sets processed at a time;
343df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
344d5b485abSSatish Balay 
345d5b485abSSatish Balay      Output:
3460e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
347d5b485abSSatish Balay                to each index set;
348d5b485abSSatish Balay       data   - pointer to the solutions
349d5b485abSSatish Balay */
350b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data)
351d5b485abSSatish Balay {
352df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
353d5b485abSSatish Balay   Mat         A  = c->A,B = c->B;
354df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
355b24ad042SBarry Smith   PetscInt    start,end,val,max,rstart,cstart,*ai,*aj;
356b24ad042SBarry Smith   PetscInt    *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
357f1af5d2fSBarry Smith   PetscBT     table_i;
358d5b485abSSatish Balay 
3593a40ed3dSBarry Smith   PetscFunctionBegin;
360899cda47SBarry Smith   rstart = c->rstartbs;
361899cda47SBarry Smith   cstart = c->cstartbs;
362d5b485abSSatish Balay   ai     = a->i;
363df36731bSBarry Smith   aj     = a->j;
364d5b485abSSatish Balay   bi     = b->i;
365df36731bSBarry Smith   bj     = b->j;
366d5b485abSSatish Balay   garray = c->garray;
367d5b485abSSatish Balay 
368d5b485abSSatish Balay 
369d5b485abSSatish Balay   for (i=0; i<imax; i++) {
370d5b485abSSatish Balay     data_i  = data[i];
371d5b485abSSatish Balay     table_i = table[i];
372d5b485abSSatish Balay     isz_i   = isz[i];
373d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
374d5b485abSSatish Balay       row   = data_i[j] - rstart;
375d5b485abSSatish Balay       start = ai[row];
376d5b485abSSatish Balay       end   = ai[row+1];
377d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
378df36731bSBarry Smith         val = aj[k] + cstart;
37926fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val;
380d5b485abSSatish Balay       }
381d5b485abSSatish Balay       start = bi[row];
382d5b485abSSatish Balay       end   = bi[row+1];
383d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
384df36731bSBarry Smith         val = garray[bj[k]];
38526fbe8dcSKarl Rupp         if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val;
386d5b485abSSatish Balay       }
387d5b485abSSatish Balay     }
388d5b485abSSatish Balay     isz[i] = isz_i;
389d5b485abSSatish Balay   }
3903a40ed3dSBarry Smith   PetscFunctionReturn(0);
391d5b485abSSatish Balay }
3924a2ae208SSatish Balay #undef __FUNCT__
3934a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive"
394d5b485abSSatish Balay /*
395df36731bSBarry Smith       MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages,
396d5b485abSSatish Balay          and return the output
397d5b485abSSatish Balay 
398d5b485abSSatish Balay          Input:
399d5b485abSSatish Balay            C    - the matrix
400d5b485abSSatish Balay            nrqr - no of messages being processed.
401d5b485abSSatish Balay            rbuf - an array of pointers to the recieved requests
402d5b485abSSatish Balay 
403d5b485abSSatish Balay          Output:
404d5b485abSSatish Balay            xdata - array of messages to be sent back
405d5b485abSSatish Balay            isz1  - size of each message
406d5b485abSSatish Balay 
407a8c7a070SBarry Smith   For better efficiency perhaps we should malloc separately each xdata[i],
408d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
4090e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of
410d5b485abSSatish Balay memory is used.
411d5b485abSSatish Balay 
412d5b485abSSatish Balay */
413b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
414d5b485abSSatish Balay {
415df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
416d5b485abSSatish Balay   Mat            A  = c->A,B = c->B;
417df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
4186849ba73SBarry Smith   PetscErrorCode ierr;
419b24ad042SBarry Smith   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
420b24ad042SBarry Smith   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
421d2a212eaSBarry Smith   PetscInt       val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
422b24ad042SBarry Smith   PetscInt       *rbuf_i,kmax,rbuf_0;
423f1af5d2fSBarry Smith   PetscBT        xtable;
424d5b485abSSatish Balay 
4253a40ed3dSBarry Smith   PetscFunctionBegin;
426df36731bSBarry Smith   Mbs    = c->Mbs;
427899cda47SBarry Smith   rstart = c->rstartbs;
428899cda47SBarry Smith   cstart = c->cstartbs;
429d5b485abSSatish Balay   ai     = a->i;
430df36731bSBarry Smith   aj     = a->j;
431d5b485abSSatish Balay   bi     = b->i;
432df36731bSBarry Smith   bj     = b->j;
433d5b485abSSatish Balay   garray = c->garray;
434d5b485abSSatish Balay 
435d5b485abSSatish Balay 
436d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
437d5b485abSSatish Balay     rbuf_i =  rbuf[i];
438d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
439d5b485abSSatish Balay     ct    += rbuf_0;
44026fbe8dcSKarl Rupp     for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j];
441d5b485abSSatish Balay   }
442d5b485abSSatish Balay 
443701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
444701b8814SBarry Smith   else        max1 = 1;
445d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
446785e854fSJed Brown   ierr         = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr);
447d5b485abSSatish Balay   ++no_malloc;
44853b8de81SBarry Smith   ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr);
449b24ad042SBarry Smith   ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr);
450d5b485abSSatish Balay 
451d5b485abSSatish Balay   ct3 = 0;
452d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
453d5b485abSSatish Balay     rbuf_i =  rbuf[i];
454d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
455d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
456d5b485abSSatish Balay     ct2    =  ct1;
457d5b485abSSatish Balay     ct3   += ct1;
458d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
4596831982aSBarry Smith       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
460d5b485abSSatish Balay       oct2 = ct2;
461d5b485abSSatish Balay       kmax = rbuf_i[2*j];
462d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
463d5b485abSSatish Balay         row = rbuf_i[ct1];
464f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
465d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
466b24ad042SBarry Smith             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
467785e854fSJed Brown             ierr         = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr);
468b24ad042SBarry Smith             ierr         = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
469606d414cSSatish Balay             ierr         = PetscFree(xdata[0]);CHKERRQ(ierr);
470d5b485abSSatish Balay             xdata[0]     = tmp;
471d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
47226fbe8dcSKarl Rupp             for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
473d5b485abSSatish Balay           }
474d5b485abSSatish Balay           xdata[i][ct2++] = row;
475d5b485abSSatish Balay           ct3++;
476d5b485abSSatish Balay         }
477d5b485abSSatish Balay       }
478d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
479d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
480d5b485abSSatish Balay         start = ai[row];
481d5b485abSSatish Balay         end   = ai[row+1];
482d5b485abSSatish Balay         for (l=start; l<end; l++) {
483df36731bSBarry Smith           val = aj[l] + cstart;
484f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
485d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
486b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
487785e854fSJed Brown               ierr         = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr);
488b24ad042SBarry Smith               ierr         = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
489606d414cSSatish Balay               ierr         = PetscFree(xdata[0]);CHKERRQ(ierr);
490d5b485abSSatish Balay               xdata[0]     = tmp;
491d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
49226fbe8dcSKarl Rupp               for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
493d5b485abSSatish Balay             }
494d5b485abSSatish Balay             xdata[i][ct2++] = val;
495d5b485abSSatish Balay             ct3++;
496d5b485abSSatish Balay           }
497d5b485abSSatish Balay         }
498d5b485abSSatish Balay         start = bi[row];
499d5b485abSSatish Balay         end   = bi[row+1];
500d5b485abSSatish Balay         for (l=start; l<end; l++) {
501df36731bSBarry Smith           val = garray[bj[l]];
502f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
503d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
504b24ad042SBarry Smith               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
505785e854fSJed Brown               ierr         = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr);
506b24ad042SBarry Smith               ierr         = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr);
507606d414cSSatish Balay               ierr         = PetscFree(xdata[0]);CHKERRQ(ierr);
508d5b485abSSatish Balay               xdata[0]     = tmp;
509d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
51026fbe8dcSKarl Rupp               for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
511d5b485abSSatish Balay             }
512d5b485abSSatish Balay             xdata[i][ct2++] = val;
513d5b485abSSatish Balay             ct3++;
514d5b485abSSatish Balay           }
515d5b485abSSatish Balay         }
516d5b485abSSatish Balay       }
517d5b485abSSatish Balay       /* Update the header*/
518d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
519d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
520d5b485abSSatish Balay     }
521d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
522d5b485abSSatish Balay     xdata[i+1]  = xdata[i] + ct2;
523d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
524d5b485abSSatish Balay   }
52594bacf5dSBarry Smith   ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr);
5261e2582c4SBarry Smith   ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr);
5273a40ed3dSBarry Smith   PetscFunctionReturn(0);
528d5b485abSSatish Balay }
529d5b485abSSatish Balay 
5304a2ae208SSatish Balay #undef __FUNCT__
5314a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ"
532b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
533d5b485abSSatish Balay {
53436f4e84dSSatish Balay   IS             *isrow_new,*iscol_new;
535cf4f527aSSatish Balay   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
5366849ba73SBarry Smith   PetscErrorCode ierr;
5374da72fa9SHong Zhang   PetscInt       nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs;
5384da72fa9SHong Zhang   PetscBool      colflag,*allcolumns,*allrows;
539a2feddc5SSatish Balay 
5403a40ed3dSBarry Smith   PetscFunctionBegin;
54129dcf524SDmitry Karpeev   /* Currently, unsorted column indices will result in inverted column indices in the resulting submatrices. */
54229dcf524SDmitry Karpeev   for (i = 0; i < ismax; ++i) {
54329dcf524SDmitry Karpeev     PetscBool sorted;
54429dcf524SDmitry Karpeev     ierr = ISSorted(iscol[i], &sorted);CHKERRQ(ierr);
545c7e6e2c7SJed Brown     if (!sorted) SETERRQ1(((PetscObject)iscol[i])->comm, PETSC_ERR_SUP, "Column index set %D not sorted", i);
54629dcf524SDmitry Karpeev   }
54729dcf524SDmitry Karpeev   /* The compression and expansion should be avoided. Doesn't point
54829dcf524SDmitry Karpeev      out errors, might change the indices, hence buggey */
549dcca6d9dSJed Brown   ierr = PetscMalloc2(ismax+1,&isrow_new,ismax+1,&iscol_new);CHKERRQ(ierr);
55005d8c843SHong Zhang   ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_new);CHKERRQ(ierr);
55105d8c843SHong Zhang   ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_new);CHKERRQ(ierr);
552cf4f527aSSatish Balay 
553307b7a18SHong Zhang   /* Check for special case: each processor gets entire matrix columns */
554dcca6d9dSJed Brown   ierr = PetscMalloc2(ismax+1,&allcolumns,ismax+1,&allrows);CHKERRQ(ierr);
555307b7a18SHong Zhang   for (i=0; i<ismax; i++) {
556307b7a18SHong Zhang     ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr);
557307b7a18SHong Zhang     ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr);
55826fbe8dcSKarl Rupp     if (colflag && ncol == C->cmap->N) allcolumns[i] = PETSC_TRUE;
55926fbe8dcSKarl Rupp     else allcolumns[i] = PETSC_FALSE;
5604da72fa9SHong Zhang 
5614da72fa9SHong Zhang     ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr);
5624da72fa9SHong Zhang     ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr);
56326fbe8dcSKarl Rupp     if (colflag && nrow == C->rmap->N) allrows[i] = PETSC_TRUE;
56426fbe8dcSKarl Rupp     else allrows[i] = PETSC_FALSE;
565307b7a18SHong Zhang   }
566307b7a18SHong Zhang 
567cf4f527aSSatish Balay   /* Allocate memory to hold all the submatrices */
568cf4f527aSSatish Balay   if (scall != MAT_REUSE_MATRIX) {
569785e854fSJed Brown     ierr = PetscMalloc1((ismax+1),submat);CHKERRQ(ierr);
570cf4f527aSSatish Balay   }
571cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
572b24ad042SBarry Smith   nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt));
573329f5518SBarry Smith   if (!nmax) nmax = 1;
574cf4f527aSSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0);
575cf4f527aSSatish Balay 
576653e4784SBarry Smith   /* Make sure every processor loops through the nstages */
577ce94432eSBarry Smith   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr);
578cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
579cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
580cf4f527aSSatish Balay     else if (pos == ismax) max_no = 0;
581cf4f527aSSatish Balay     else                   max_no = ismax-pos;
5824da72fa9SHong Zhang     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr);
583cf4f527aSSatish Balay     pos += max_no;
584cf4f527aSSatish Balay   }
58536f4e84dSSatish Balay 
58636f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
5876bf464f9SBarry Smith     ierr = ISDestroy(&isrow_new[i]);CHKERRQ(ierr);
5886bf464f9SBarry Smith     ierr = ISDestroy(&iscol_new[i]);CHKERRQ(ierr);
58936f4e84dSSatish Balay   }
59023bdbc58SBarry Smith   ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr);
5914da72fa9SHong Zhang   ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr);
5923a40ed3dSBarry Smith   PetscFunctionReturn(0);
593a2feddc5SSatish Balay }
594a2feddc5SSatish Balay 
595233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
5964a2ae208SSatish Balay #undef __FUNCT__
5974a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc"
598e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank)
599233c2ff6SSatish Balay {
600e005ede5SBarry Smith   PetscInt       nGlobalNd = proc_gnode[size];
6014dc2109aSBarry Smith   PetscMPIInt    fproc;
6024dc2109aSBarry Smith   PetscErrorCode ierr;
603233c2ff6SSatish Balay 
604233c2ff6SSatish Balay   PetscFunctionBegin;
6054dc2109aSBarry Smith   ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr);
60623ce1328SBarry Smith   if (fproc > size) fproc = size;
607e005ede5SBarry Smith   while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) {
608e005ede5SBarry Smith     if (row < proc_gnode[fproc]) fproc--;
609233c2ff6SSatish Balay     else                         fproc++;
610233c2ff6SSatish Balay   }
611e005ede5SBarry Smith   *rank = fproc;
612233c2ff6SSatish Balay   PetscFunctionReturn(0);
613233c2ff6SSatish Balay }
614233c2ff6SSatish Balay #endif
615233c2ff6SSatish Balay 
616a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
617b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */
6184a2ae208SSatish Balay #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local"
6204da72fa9SHong Zhang PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats)
621a2feddc5SSatish Balay {
622df36731bSBarry Smith   Mat_MPIBAIJ    *c = (Mat_MPIBAIJ*)C->data;
623cf4f527aSSatish Balay   Mat            A  = c->A;
624df36731bSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
6255d0c19d7SBarry Smith   const PetscInt **irow,**icol,*irow_i;
6265d0c19d7SBarry Smith   PetscInt       *nrow,*ncol,*w3,*w4,start;
6276849ba73SBarry Smith   PetscErrorCode ierr;
6289405f653SBarry Smith   PetscMPIInt    size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc;
6299405f653SBarry Smith   PetscInt       **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row;
630b24ad042SBarry Smith   PetscInt       nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol;
631b24ad042SBarry Smith   PetscInt       **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
6325d0c19d7SBarry Smith   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
633052f0c41SBarry Smith   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
634d0f46423SBarry Smith   PetscInt       bs     =C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
635899cda47SBarry Smith   PetscInt       cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
636899cda47SBarry Smith   PetscInt       *bmap  = c->garray,ctmp,rstart=c->rstartbs;
6377a868f3eSHong Zhang   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3;
6387a868f3eSHong Zhang   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3;
639d5b485abSSatish Balay   MPI_Comm       comm;
640ace3abfcSBarry Smith   PetscBool      flag;
641b24ad042SBarry Smith   PetscMPIInt    *onodes1,*olengths1;
6427a868f3eSHong Zhang   PetscBool      ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */
64326fbe8dcSKarl Rupp 
6447a868f3eSHong Zhang   /* variables below are used for the matrix numerical values - case of !ijonly */
6457a868f3eSHong Zhang   MPI_Request *r_waits4,*s_waits4;
6467a868f3eSHong Zhang   MPI_Status  *r_status4,*s_status4;
6470298fd71SBarry Smith   MatScalar   **rbuf4,**sbuf_aa,*vals,*mat_a = NULL,*sbuf_aa_i,*vworkA = NULL,*vworkB = NULL;
6487a868f3eSHong Zhang   MatScalar   *a_a=a->a,*b_a=b->a;
649c7dd2924SSatish Balay 
65080d608c0SSatish Balay #if defined(PETSC_USE_CTABLE)
651b24ad042SBarry Smith   PetscInt   tt;
6520298fd71SBarry Smith   PetscTable *rmap,*cmap,rmap_i,cmap_i=NULL;
653233c2ff6SSatish Balay #else
6540298fd71SBarry Smith   PetscInt **cmap,*cmap_i=NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
655233c2ff6SSatish Balay #endif
656d5b485abSSatish Balay 
6573a40ed3dSBarry Smith   PetscFunctionBegin;
658ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
6597adad957SLisandro Dalcin   tag0 = ((PetscObject)C)->tag;
660d5b485abSSatish Balay   size = c->size;
661d5b485abSSatish Balay   rank = c->rank;
662d5b485abSSatish Balay 
66334e6ae18SSatish Balay   /* Get some new tags to keep the communication clean */
66434e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
66534e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
66634e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
66734e6ae18SSatish Balay 
668052f0c41SBarry Smith #if defined(PETSC_USE_CTABLE)
669dcca6d9dSJed Brown   ierr = PetscMalloc4(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol);CHKERRQ(ierr);
670052f0c41SBarry Smith #else
671dcca6d9dSJed Brown   ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,Mbs+1,&rtable);CHKERRQ(ierr);
672d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
673d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
67482a7e548SBarry Smith     jmax = C->rmap->range[i+1]/bs;
67526fbe8dcSKarl Rupp     for (; j<jmax; j++) rtable[j] = i;
676d5b485abSSatish Balay   }
677233c2ff6SSatish Balay #endif
678233c2ff6SSatish Balay 
679233c2ff6SSatish Balay   for (i=0; i<ismax; i++) {
6804da72fa9SHong Zhang     if (allrows[i]) {
6810298fd71SBarry Smith       irow[i] = NULL;
6824da72fa9SHong Zhang       nrow[i] = C->rmap->N/bs;
6834da72fa9SHong Zhang     } else {
684233c2ff6SSatish Balay       ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
685233c2ff6SSatish Balay       ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
6864da72fa9SHong Zhang     }
6874da72fa9SHong Zhang 
688307b7a18SHong Zhang     if (allcolumns[i]) {
6890298fd71SBarry Smith       icol[i] = NULL;
690307b7a18SHong Zhang       ncol[i] = C->cmap->N/bs;
691307b7a18SHong Zhang     } else {
692307b7a18SHong Zhang       ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
693233c2ff6SSatish Balay       ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
694233c2ff6SSatish Balay     }
695307b7a18SHong Zhang   }
696d5b485abSSatish Balay 
697d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg,and buffer space
698d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
699dcca6d9dSJed Brown   ierr = PetscMalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr);
70023bdbc58SBarry Smith   ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
70123bdbc58SBarry Smith   ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
70223bdbc58SBarry Smith   ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr);
703d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
704b24ad042SBarry Smith     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/
705d5b485abSSatish Balay     jmax   = nrow[i];
706d5b485abSSatish Balay     irow_i = irow[i];
707d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
70826fbe8dcSKarl Rupp       if (allrows[i]) row = j;
70926fbe8dcSKarl Rupp       else row = irow_i[j];
71026fbe8dcSKarl Rupp 
711233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
712899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
713233c2ff6SSatish Balay #else
714d5b485abSSatish Balay       proc = rtable[row];
715233c2ff6SSatish Balay #endif
716d5b485abSSatish Balay       w4[proc]++;
717d5b485abSSatish Balay     }
718d5b485abSSatish Balay     for (j=0; j<size; j++) {
719d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
720d5b485abSSatish Balay     }
721d5b485abSSatish Balay   }
722d5b485abSSatish Balay 
723d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
724e757655aSSatish Balay   msz      = 0;              /* total mesg length for all proc */
725d5b485abSSatish Balay   w1[rank] = 0;              /* no mesg sent to intself */
726d5b485abSSatish Balay   w3[rank] = 0;
727d5b485abSSatish Balay   for (i=0; i<size; i++) {
728d5b485abSSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
729d5b485abSSatish Balay   }
730785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&pa);CHKERRQ(ierr); /*(proc -array)*/
731d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
732d5b485abSSatish Balay     if (w1[i]) { pa[j] = i; j++; }
733d5b485abSSatish Balay   }
734d5b485abSSatish Balay 
735d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
736d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
737d5b485abSSatish Balay     j     = pa[i];
738d5b485abSSatish Balay     w1[j] += w2[j] + 2* w3[j];
739d5b485abSSatish Balay     msz   += w1[j];
740d5b485abSSatish Balay   }
741d5b485abSSatish Balay 
742c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
743a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
744c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
745d5b485abSSatish Balay 
746c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
747c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
748c7dd2924SSatish Balay 
749c7dd2924SSatish Balay   ierr = PetscFree(onodes1);CHKERRQ(ierr);
750c7dd2924SSatish Balay   ierr = PetscFree(olengths1);CHKERRQ(ierr);
751d5b485abSSatish Balay 
752d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
753dcca6d9dSJed Brown   ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr);
75423bdbc58SBarry Smith   ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr);
75523bdbc58SBarry Smith   ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr);
756d5b485abSSatish Balay   {
757b24ad042SBarry Smith     PetscInt *iptr = tmp,ict = 0;
758d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
759d5b485abSSatish Balay       j        = pa[i];
760d5b485abSSatish Balay       iptr    += ict;
761d5b485abSSatish Balay       sbuf1[j] = iptr;
762d5b485abSSatish Balay       ict      = w1[j];
763d5b485abSSatish Balay     }
764d5b485abSSatish Balay   }
765d5b485abSSatish Balay 
766d5b485abSSatish Balay   /* Form the outgoing messages */
767d5b485abSSatish Balay   /* Initialise the header space */
768d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
769d5b485abSSatish Balay     j           = pa[i];
770d5b485abSSatish Balay     sbuf1[j][0] = 0;
771b24ad042SBarry Smith     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
772d5b485abSSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
773d5b485abSSatish Balay   }
774d5b485abSSatish Balay 
775d5b485abSSatish Balay   /* Parse the isrow and copy data into outbuf */
776d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
777b24ad042SBarry Smith     ierr   = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
778d5b485abSSatish Balay     irow_i = irow[i];
779d5b485abSSatish Balay     jmax   = nrow[i];
780d5b485abSSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
78126fbe8dcSKarl Rupp       if (allrows[i]) row = j;
78226fbe8dcSKarl Rupp       else row = irow_i[j];
78326fbe8dcSKarl Rupp 
784233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
785899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
786233c2ff6SSatish Balay #else
787d5b485abSSatish Balay       proc = rtable[row];
788233c2ff6SSatish Balay #endif
789d5b485abSSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
790d5b485abSSatish Balay         ctr[proc]++;
791d5b485abSSatish Balay         *ptr[proc] = row;
792d5b485abSSatish Balay         ptr[proc]++;
793d5b485abSSatish Balay       }
794d5b485abSSatish Balay     }
795d5b485abSSatish Balay     /* Update the headers for the current IS */
796d5b485abSSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
79706ef35abSSatish Balay       if ((ctr_j = ctr[j])) {
798d5b485abSSatish Balay         sbuf1_j        = sbuf1[j];
799d5b485abSSatish Balay         k              = ++sbuf1_j[0];
800d5b485abSSatish Balay         sbuf1_j[2*k]   = ctr_j;
801d5b485abSSatish Balay         sbuf1_j[2*k-1] = i;
802d5b485abSSatish Balay       }
803d5b485abSSatish Balay     }
804d5b485abSSatish Balay   }
805d5b485abSSatish Balay 
806d5b485abSSatish Balay   /*  Now  post the sends */
807785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&s_waits1);CHKERRQ(ierr);
808d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
809d5b485abSSatish Balay     j    = pa[i];
810b24ad042SBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
811d5b485abSSatish Balay   }
812d5b485abSSatish Balay 
813d5b485abSSatish Balay   /* Post Recieves to capture the buffer size */
814785e854fSJed Brown   ierr     = PetscMalloc1((nrqs+1),&r_waits2);CHKERRQ(ierr);
815785e854fSJed Brown   ierr     = PetscMalloc1((nrqs+1),&rbuf2);CHKERRQ(ierr);
816d5b485abSSatish Balay   rbuf2[0] = tmp + msz;
817d5b485abSSatish Balay   for (i=1; i<nrqs; ++i) {
818d5b485abSSatish Balay     j        = pa[i];
819d5b485abSSatish Balay     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
820d5b485abSSatish Balay   }
821d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
822d5b485abSSatish Balay     j    = pa[i];
823b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
824d5b485abSSatish Balay   }
825d5b485abSSatish Balay 
826d5b485abSSatish Balay   /* Send to other procs the buf size they should allocate */
827d5b485abSSatish Balay 
828d5b485abSSatish Balay   /* Receive messages*/
829785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&s_waits2);CHKERRQ(ierr);
830785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&r_status1);CHKERRQ(ierr);
831dcca6d9dSJed Brown   ierr = PetscMalloc3(nrqr+1,&sbuf2,nrqr,&req_size,nrqr,&req_source);CHKERRQ(ierr);
832d5b485abSSatish Balay   {
833df36731bSBarry Smith     Mat_SeqBAIJ *sA  = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
834b24ad042SBarry Smith     PetscInt    *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
835d5b485abSSatish Balay 
836d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
837999d9058SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
83826fbe8dcSKarl Rupp 
839999d9058SBarry Smith       req_size[idex] = 0;
840999d9058SBarry Smith       rbuf1_i        = rbuf1[idex];
841d5b485abSSatish Balay       start          = 2*rbuf1_i[0] + 1;
842b24ad042SBarry Smith       ierr           = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
843785e854fSJed Brown       ierr           = PetscMalloc1(end,&sbuf2[idex]);CHKERRQ(ierr);
844999d9058SBarry Smith       sbuf2_i        = sbuf2[idex];
845d5b485abSSatish Balay       for (j=start; j<end; j++) {
846d5b485abSSatish Balay         id              = rbuf1_i[j] - rstart;
847d5b485abSSatish Balay         ncols           = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
848d5b485abSSatish Balay         sbuf2_i[j]      = ncols;
849999d9058SBarry Smith         req_size[idex] += ncols;
850d5b485abSSatish Balay       }
851999d9058SBarry Smith       req_source[idex] = r_status1[i].MPI_SOURCE;
852d5b485abSSatish Balay       /* form the header */
853999d9058SBarry Smith       sbuf2_i[0] = req_size[idex];
85426fbe8dcSKarl Rupp       for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
855b24ad042SBarry Smith       ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
856d5b485abSSatish Balay     }
857d5b485abSSatish Balay   }
858606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
859606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
860d5b485abSSatish Balay 
861d5b485abSSatish Balay   /*  recv buffer sizes */
862d5b485abSSatish Balay   /* Receive messages*/
863785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&rbuf3);CHKERRQ(ierr);
864785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&r_waits3);CHKERRQ(ierr);
865785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&r_status2);CHKERRQ(ierr);
8667a868f3eSHong Zhang   if (!ijonly) {
867785e854fSJed Brown     ierr = PetscMalloc1((nrqs+1),&rbuf4);CHKERRQ(ierr);
868785e854fSJed Brown     ierr = PetscMalloc1((nrqs+1),&r_waits4);CHKERRQ(ierr);
8697a868f3eSHong Zhang   }
870d5b485abSSatish Balay 
871d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
872999d9058SBarry Smith     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
873785e854fSJed Brown     ierr = PetscMalloc1(rbuf2[idex][0],&rbuf3[idex]);CHKERRQ(ierr);
874b24ad042SBarry Smith     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
8757a868f3eSHong Zhang     if (!ijonly) {
876785e854fSJed Brown       ierr = PetscMalloc1(rbuf2[idex][0]*bs2,&rbuf4[idex]);CHKERRQ(ierr);
877b24ad042SBarry Smith       ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
878d5b485abSSatish Balay     }
8797a868f3eSHong Zhang   }
880606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
881606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
882d5b485abSSatish Balay 
883d5b485abSSatish Balay   /* Wait on sends1 and sends2 */
884785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&s_status1);CHKERRQ(ierr);
885785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&s_status2);CHKERRQ(ierr);
886d5b485abSSatish Balay 
8870c468ba9SBarry Smith   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);}
8880c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);}
889606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
890606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
891606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
892606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
893d5b485abSSatish Balay 
894d5b485abSSatish Balay   /* Now allocate buffers for a->j, and send them off */
895785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&sbuf_aj);CHKERRQ(ierr);
896d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
897785e854fSJed Brown   ierr = PetscMalloc1((j+1),&sbuf_aj[0]);CHKERRQ(ierr);
898d5b485abSSatish Balay   for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
899d5b485abSSatish Balay 
900785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&s_waits3);CHKERRQ(ierr);
901d5b485abSSatish Balay   {
902d5b485abSSatish Balay     for (i=0; i<nrqr; i++) {
903d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
904d5b485abSSatish Balay       sbuf_aj_i = sbuf_aj[i];
905d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0] + 1;
906d5b485abSSatish Balay       ct2       = 0;
907d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
908d5b485abSSatish Balay         kmax = rbuf1[i][2*j];
909d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
910e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
911d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
912d5b485abSSatish Balay           ncols  = nzA + nzB;
913d5b485abSSatish Balay           cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
914d5b485abSSatish Balay 
915d5b485abSSatish Balay           /* load the column indices for this row into cols*/
916d5b485abSSatish Balay           cols = sbuf_aj_i + ct2;
917d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
918d5b485abSSatish Balay             if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp;
919d5b485abSSatish Balay             else break;
920d5b485abSSatish Balay           }
921d5b485abSSatish Balay           imark = l;
922d5b485abSSatish Balay           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
923d5b485abSSatish Balay           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
924d5b485abSSatish Balay           ct2 += ncols;
925d5b485abSSatish Balay         }
926d5b485abSSatish Balay       }
927b24ad042SBarry Smith       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
928d5b485abSSatish Balay     }
929d5b485abSSatish Balay   }
930785e854fSJed Brown   ierr = PetscMalloc1((nrqs+1),&r_status3);CHKERRQ(ierr);
931785e854fSJed Brown   ierr = PetscMalloc1((nrqr+1),&s_status3);CHKERRQ(ierr);
932d5b485abSSatish Balay 
933d5b485abSSatish Balay   /* Allocate buffers for a->a, and send them off */
9347a868f3eSHong Zhang   if (!ijonly) {
935785e854fSJed Brown     ierr = PetscMalloc1((nrqr+1),&sbuf_aa);CHKERRQ(ierr);
936d5b485abSSatish Balay     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
937785e854fSJed Brown     ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr);
938a2feddc5SSatish Balay     for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
939d5b485abSSatish Balay 
940785e854fSJed Brown     ierr = PetscMalloc1((nrqr+1),&s_waits4);CHKERRQ(ierr);
941d5b485abSSatish Balay     {
942d5b485abSSatish Balay       for (i=0; i<nrqr; i++) {
943d5b485abSSatish Balay         rbuf1_i   = rbuf1[i];
944d5b485abSSatish Balay         sbuf_aa_i = sbuf_aa[i];
945d5b485abSSatish Balay         ct1       = 2*rbuf1_i[0]+1;
946d5b485abSSatish Balay         ct2       = 0;
947d5b485abSSatish Balay         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
948d5b485abSSatish Balay           kmax = rbuf1_i[2*j];
949d5b485abSSatish Balay           for (k=0; k<kmax; k++,ct1++) {
950e8e0db45SSatish Balay             row    = rbuf1_i[ct1] - rstart;
951d5b485abSSatish Balay             nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
952d5b485abSSatish Balay             ncols  = nzA + nzB;
953d5b485abSSatish Balay             cworkA = a_j + a_i[row];     cworkB = b_j + b_i[row];
954a2feddc5SSatish Balay             vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2;
955d5b485abSSatish Balay 
956d5b485abSSatish Balay             /* load the column values for this row into vals*/
9575b83ace0SSatish Balay             vals = sbuf_aa_i+ct2*bs2;
958d5b485abSSatish Balay             for (l=0; l<nzB; l++) {
9593a40ed3dSBarry Smith               if ((bmap[cworkB[l]]) < cstart) {
9603eda8832SBarry Smith                 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
96126fbe8dcSKarl Rupp               } else break;
962d5b485abSSatish Balay             }
963d5b485abSSatish Balay             imark = l;
9643a40ed3dSBarry Smith             for (l=0; l<nzA; l++) {
9653eda8832SBarry Smith               ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9663a40ed3dSBarry Smith             }
9673a40ed3dSBarry Smith             for (l=imark; l<nzB; l++) {
9683eda8832SBarry Smith               ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9693a40ed3dSBarry Smith             }
970d5b485abSSatish Balay             ct2 += ncols;
971d5b485abSSatish Balay           }
972d5b485abSSatish Balay         }
9733eda8832SBarry Smith         ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
974d5b485abSSatish Balay       }
975d5b485abSSatish Balay     }
976785e854fSJed Brown     ierr = PetscMalloc1((nrqs+1),&r_status4);CHKERRQ(ierr);
977785e854fSJed Brown     ierr = PetscMalloc1((nrqr+1),&s_status4);CHKERRQ(ierr);
9787a868f3eSHong Zhang   }
979533163c2SBarry Smith   ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr);
980606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
981d5b485abSSatish Balay 
982d5b485abSSatish Balay   /* Form the matrix */
983307b7a18SHong Zhang   /* create col map: global col of C -> local col of submatrices */
984d5b485abSSatish Balay   {
9855d0c19d7SBarry Smith     const PetscInt *icol_i;
986233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
987785e854fSJed Brown     ierr = PetscMalloc1((1+ismax),&cmap);CHKERRQ(ierr);
988ff0794f7SSatish Balay     for (i=0; i<ismax; i++) {
989307b7a18SHong Zhang       if (!allcolumns[i]) {
9908fa52d88SHong Zhang         ierr   = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr);
991307b7a18SHong Zhang         jmax   = ncol[i];
992307b7a18SHong Zhang         icol_i = icol[i];
9938fa52d88SHong Zhang         cmap_i = cmap[i];
994307b7a18SHong Zhang         for (j=0; j<jmax; j++) {
9953861aac3SJed Brown           ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
996307b7a18SHong Zhang         }
997307b7a18SHong Zhang       } else {
9980298fd71SBarry Smith         cmap[i] = NULL;
999307b7a18SHong Zhang       }
1000233c2ff6SSatish Balay     }
1001233c2ff6SSatish Balay #else
1002785e854fSJed Brown     ierr = PetscMalloc1(ismax,&cmap);CHKERRQ(ierr);
1003d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
10048fa52d88SHong Zhang       if (!allcolumns[i]) {
1005785e854fSJed Brown         ierr   = PetscMalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr);
10068fa52d88SHong Zhang         ierr   = PetscMemzero(cmap[i],c->Nbs*sizeof(PetscInt));CHKERRQ(ierr);
1007d5b485abSSatish Balay         jmax   = ncol[i];
1008d5b485abSSatish Balay         icol_i = icol[i];
1009d5b485abSSatish Balay         cmap_i = cmap[i];
1010d5b485abSSatish Balay         for (j=0; j<jmax; j++) {
1011d5b485abSSatish Balay           cmap_i[icol_i[j]] = j+1;
1012d5b485abSSatish Balay         }
10138fa52d88SHong Zhang       } else { /* allcolumns[i] */
10140298fd71SBarry Smith         cmap[i] = NULL;
10158fa52d88SHong Zhang       }
1016d5b485abSSatish Balay     }
1017307b7a18SHong Zhang #endif
1018d5b485abSSatish Balay   }
1019d5b485abSSatish Balay 
1020d5b485abSSatish Balay   /* Create lens which is required for MatCreate... */
102126fbe8dcSKarl Rupp   for (i=0,j=0; i<ismax; i++) j += nrow[i];
1022785e854fSJed Brown   ierr    = PetscMalloc1((1+ismax)*sizeof(PetscInt*)+ j,&lens);CHKERRQ(ierr);
1023b24ad042SBarry Smith   lens[0] = (PetscInt*)(lens + ismax);
1024b24ad042SBarry Smith   ierr    = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr);
102526fbe8dcSKarl Rupp   for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];
1026d5b485abSSatish Balay 
1027d5b485abSSatish Balay   /* Update lens from local data */
1028d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1029d5b485abSSatish Balay     jmax = nrow[i];
10308fa52d88SHong Zhang     if (!allcolumns[i]) cmap_i = cmap[i];
1031d5b485abSSatish Balay     irow_i = irow[i];
1032d5b485abSSatish Balay     lens_i = lens[i];
1033d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
103426fbe8dcSKarl Rupp       if (allrows[i]) row = j;
103526fbe8dcSKarl Rupp       else row = irow_i[j];
103626fbe8dcSKarl Rupp 
1037233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
1038899cda47SBarry Smith       ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1039233c2ff6SSatish Balay #else
1040d5b485abSSatish Balay       proc = rtable[row];
1041233c2ff6SSatish Balay #endif
1042d5b485abSSatish Balay       if (proc == rank) {
104336f4e84dSSatish Balay         /* Get indices from matA and then from matB */
104436f4e84dSSatish Balay         row    = row - rstart;
104536f4e84dSSatish Balay         nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
104636f4e84dSSatish Balay         cworkA =  a_j + a_i[row]; cworkB = b_j + b_i[row];
1047307b7a18SHong Zhang         if (!allcolumns[i]) {
1048233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
1049233c2ff6SSatish Balay           for (k=0; k<nzA; k++) {
10508fa52d88SHong Zhang             ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
105126fbe8dcSKarl Rupp             if (tt) lens_i[j]++;
1052233c2ff6SSatish Balay           }
1053233c2ff6SSatish Balay           for (k=0; k<nzB; k++) {
10548fa52d88SHong Zhang             ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
105526fbe8dcSKarl Rupp             if (tt) lens_i[j]++;
1056233c2ff6SSatish Balay           }
1057307b7a18SHong Zhang 
1058233c2ff6SSatish Balay #else
1059ca161407SBarry Smith           for (k=0; k<nzA; k++) {
106026fbe8dcSKarl Rupp             if (cmap_i[cstart + cworkA[k]]) lens_i[j]++;
1061ca161407SBarry Smith           }
1062ca161407SBarry Smith           for (k=0; k<nzB; k++) {
106326fbe8dcSKarl Rupp             if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++;
1064d5b485abSSatish Balay           }
1065233c2ff6SSatish Balay #endif
1066307b7a18SHong Zhang         } else { /* allcolumns */
1067307b7a18SHong Zhang           lens_i[j] = nzA + nzB;
1068307b7a18SHong Zhang         }
1069d5b485abSSatish Balay       }
1070a2feddc5SSatish Balay     }
1071ca161407SBarry Smith   }
1072233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
1073d5b485abSSatish Balay   /* Create row map*/
1074785e854fSJed Brown   ierr = PetscMalloc1((1+ismax),&rmap);CHKERRQ(ierr);
1075ff0794f7SSatish Balay   for (i=0; i<ismax; i++) {
10768fa52d88SHong Zhang     ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr);
1077233c2ff6SSatish Balay   }
1078233c2ff6SSatish Balay #else
1079233c2ff6SSatish Balay   /* Create row map*/
1080785e854fSJed Brown   ierr    = PetscMalloc1((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs,&rmap);CHKERRQ(ierr);
1081b24ad042SBarry Smith   rmap[0] = (PetscInt*)(rmap + ismax);
1082b24ad042SBarry Smith   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr);
108326fbe8dcSKarl Rupp   for (i=1; i<ismax; i++) rmap[i] = rmap[i-1] + Mbs;
1084233c2ff6SSatish Balay #endif
1085d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1086d5b485abSSatish Balay     irow_i = irow[i];
1087d5b485abSSatish Balay     jmax   = nrow[i];
1088233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
10898fa52d88SHong Zhang     rmap_i = rmap[i];
1090233c2ff6SSatish Balay     for (j=0; j<jmax; j++) {
10914da72fa9SHong Zhang       if (allrows[i]) {
10923861aac3SJed Brown         ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
10934da72fa9SHong Zhang       } else {
10943861aac3SJed Brown         ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr);
1095233c2ff6SSatish Balay       }
10964da72fa9SHong Zhang     }
1097233c2ff6SSatish Balay #else
1098233c2ff6SSatish Balay     rmap_i = rmap[i];
1099d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
110026fbe8dcSKarl Rupp       if (allrows[i]) rmap_i[j] = j;
110126fbe8dcSKarl Rupp       else rmap_i[irow_i[j]] = j;
11024da72fa9SHong Zhang     }
1103233c2ff6SSatish Balay #endif
1104d5b485abSSatish Balay   }
1105d5b485abSSatish Balay 
1106d5b485abSSatish Balay   /* Update lens from offproc data */
1107d5b485abSSatish Balay   {
1108b24ad042SBarry Smith     PetscInt    *rbuf2_i,*rbuf3_i,*sbuf1_i;
1109b24ad042SBarry Smith     PetscMPIInt ii;
1110d5b485abSSatish Balay 
1111d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1112b24ad042SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr);
1113b24ad042SBarry Smith       idex    = pa[ii];
1114999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1115d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1116d5b485abSSatish Balay       ct1     = 2*jmax+1;
1117d5b485abSSatish Balay       ct2     = 0;
1118b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1119b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
1120d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1121d5b485abSSatish Balay         is_no  = sbuf1_i[2*j-1];
1122d5b485abSSatish Balay         max1   = sbuf1_i[2*j];
1123d5b485abSSatish Balay         lens_i = lens[is_no];
11248fa52d88SHong Zhang         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1125d5b485abSSatish Balay         rmap_i = rmap[is_no];
1126d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1127233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
11288fa52d88SHong Zhang           ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1129233c2ff6SSatish Balay           row--;
113026fbe8dcSKarl Rupp           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1131233c2ff6SSatish Balay #else
1132d5b485abSSatish Balay           row = rmap_i[sbuf1_i[ct1]];  /* the val in the new matrix to be */
1133233c2ff6SSatish Balay #endif
1134d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1135d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1136307b7a18SHong Zhang             if (!allcolumns[is_no]) {
1137233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
11388fa52d88SHong Zhang               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
113926fbe8dcSKarl Rupp               if (tt) lens_i[row]++;
1140233c2ff6SSatish Balay #else
114126fbe8dcSKarl Rupp               if (cmap_i[rbuf3_i[ct2]]) lens_i[row]++;
1142233c2ff6SSatish Balay #endif
1143307b7a18SHong Zhang             } else { /* allcolumns */
1144307b7a18SHong Zhang               lens_i[row]++;
1145307b7a18SHong Zhang             }
1146d5b485abSSatish Balay           }
1147d5b485abSSatish Balay         }
1148d5b485abSSatish Balay       }
1149d5b485abSSatish Balay     }
1150d5b485abSSatish Balay   }
1151606d414cSSatish Balay   ierr = PetscFree(r_status3);CHKERRQ(ierr);
1152606d414cSSatish Balay   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
11530c468ba9SBarry Smith   if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);}
1154606d414cSSatish Balay   ierr = PetscFree(s_status3);CHKERRQ(ierr);
1155606d414cSSatish Balay   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1156d5b485abSSatish Balay 
1157d5b485abSSatish Balay   /* Create the submatrices */
1158d5b485abSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
11597a868f3eSHong Zhang     if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet");
1160d5b485abSSatish Balay     /*
1161d5b485abSSatish Balay         Assumes new rows are same length as the old rows, hence bug!
1162d5b485abSSatish Balay     */
1163d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1164df36731bSBarry Smith       mat = (Mat_SeqBAIJ*)(submats[i]->data);
1165e7e72b3dSBarry Smith       if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1166b24ad042SBarry Smith       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
1167e7e72b3dSBarry Smith       if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1168d5b485abSSatish Balay       /* Initial matrix as if empty */
1169b24ad042SBarry Smith       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr);
117026fbe8dcSKarl Rupp 
1171d5f3da31SBarry Smith       submats[i]->factortype = C->factortype;
1172d5b485abSSatish Balay     }
1173ca161407SBarry Smith   } else {
11747a868f3eSHong Zhang     PetscInt bs_tmp;
117526fbe8dcSKarl Rupp     if (ijonly) bs_tmp = 1;
117626fbe8dcSKarl Rupp     else        bs_tmp = bs;
1177d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1178f69a0ea3SMatthew Knepley       ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr);
11797a868f3eSHong Zhang       ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr);
11807adad957SLisandro Dalcin       ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr);
11817a868f3eSHong Zhang       ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr);
11827a868f3eSHong Zhang       ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */
1183d5b485abSSatish Balay     }
1184d5b485abSSatish Balay   }
1185d5b485abSSatish Balay 
1186d5b485abSSatish Balay   /* Assemble the matrices */
1187d5b485abSSatish Balay   /* First assemble the local rows */
1188d5b485abSSatish Balay   {
1189b24ad042SBarry Smith     PetscInt  ilen_row,*imat_ilen,*imat_j,*imat_i;
11900298fd71SBarry Smith     MatScalar *imat_a = NULL;
1191d5b485abSSatish Balay 
1192d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1193df36731bSBarry Smith       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1194d5b485abSSatish Balay       imat_ilen = mat->ilen;
1195d5b485abSSatish Balay       imat_j    = mat->j;
1196d5b485abSSatish Balay       imat_i    = mat->i;
11977a868f3eSHong Zhang       if (!ijonly) imat_a = mat->a;
11988fa52d88SHong Zhang       if (!allcolumns[i]) cmap_i = cmap[i];
1199d5b485abSSatish Balay       rmap_i = rmap[i];
1200d5b485abSSatish Balay       irow_i = irow[i];
1201d5b485abSSatish Balay       jmax   = nrow[i];
1202d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
120326fbe8dcSKarl Rupp         if (allrows[i]) row = j;
120426fbe8dcSKarl Rupp         else row = irow_i[j];
120526fbe8dcSKarl Rupp 
1206233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
1207899cda47SBarry Smith         ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr);
1208233c2ff6SSatish Balay #else
1209d5b485abSSatish Balay         proc = rtable[row];
1210233c2ff6SSatish Balay #endif
1211d5b485abSSatish Balay         if (proc == rank) {
121236f4e84dSSatish Balay           row    = row - rstart;
121336f4e84dSSatish Balay           nzA    = a_i[row+1] - a_i[row];
121436f4e84dSSatish Balay           nzB    = b_i[row+1] - b_i[row];
121536f4e84dSSatish Balay           cworkA = a_j + a_i[row];
121636f4e84dSSatish Balay           cworkB = b_j + b_i[row];
12177a868f3eSHong Zhang           if (!ijonly) {
121836f4e84dSSatish Balay             vworkA = a_a + a_i[row]*bs2;
121936f4e84dSSatish Balay             vworkB = b_a + b_i[row]*bs2;
12207a868f3eSHong Zhang           }
1221233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
12228fa52d88SHong Zhang           ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr);
1223233c2ff6SSatish Balay           row--;
122426fbe8dcSKarl Rupp           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1225233c2ff6SSatish Balay #else
122636f4e84dSSatish Balay           row = rmap_i[row + rstart];
1227233c2ff6SSatish Balay #endif
1228df36731bSBarry Smith           mat_i = imat_i[row];
12297a868f3eSHong Zhang           if (!ijonly) mat_a = imat_a + mat_i*bs2;
1230d5b485abSSatish Balay           mat_j    = imat_j + mat_i;
123136f4e84dSSatish Balay           ilen_row = imat_ilen[row];
123236f4e84dSSatish Balay 
123336f4e84dSSatish Balay           /* load the column indices for this row into cols*/
1234307b7a18SHong Zhang           if (!allcolumns[i]) {
123536f4e84dSSatish Balay             for (l=0; l<nzB; l++) {
123636f4e84dSSatish Balay               if ((ctmp = bmap[cworkB[l]]) < cstart) {
1237233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
12388fa52d88SHong Zhang                 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr);
1239233c2ff6SSatish Balay                 if (tcol) {
1240233c2ff6SSatish Balay #else
124136f4e84dSSatish Balay                 if ((tcol = cmap_i[ctmp])) {
1242233c2ff6SSatish Balay #endif
1243df36731bSBarry Smith                   *mat_j++ = tcol - 1;
12443eda8832SBarry Smith                   ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1245549d3d68SSatish Balay                   mat_a   += bs2;
1246d5b485abSSatish Balay                   ilen_row++;
1247d5b485abSSatish Balay                 }
1248ca161407SBarry Smith               } else break;
124936f4e84dSSatish Balay             }
125036f4e84dSSatish Balay             imark = l;
125136f4e84dSSatish Balay             for (l=0; l<nzA; l++) {
1252233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
12538fa52d88SHong Zhang               ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1254233c2ff6SSatish Balay               if (tcol) {
1255233c2ff6SSatish Balay #else
125636f4e84dSSatish Balay               if ((tcol = cmap_i[cstart + cworkA[l]])) {
1257233c2ff6SSatish Balay #endif
125836f4e84dSSatish Balay                 *mat_j++ = tcol - 1;
12597a868f3eSHong Zhang                 if (!ijonly) {
12603eda8832SBarry Smith                   ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1261549d3d68SSatish Balay                   mat_a += bs2;
12627a868f3eSHong Zhang                 }
126336f4e84dSSatish Balay                 ilen_row++;
126436f4e84dSSatish Balay               }
126536f4e84dSSatish Balay             }
126636f4e84dSSatish Balay             for (l=imark; l<nzB; l++) {
1267233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
12688fa52d88SHong Zhang               ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1269233c2ff6SSatish Balay               if (tcol) {
1270233c2ff6SSatish Balay #else
127136f4e84dSSatish Balay               if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1272233c2ff6SSatish Balay #endif
127336f4e84dSSatish Balay                 *mat_j++ = tcol - 1;
12747a868f3eSHong Zhang                 if (!ijonly) {
12753eda8832SBarry Smith                   ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1276549d3d68SSatish Balay                   mat_a += bs2;
12777a868f3eSHong Zhang                 }
127836f4e84dSSatish Balay                 ilen_row++;
127936f4e84dSSatish Balay               }
128036f4e84dSSatish Balay             }
1281307b7a18SHong Zhang           } else { /* allcolumns */
1282307b7a18SHong Zhang             for (l=0; l<nzB; l++) {
1283307b7a18SHong Zhang               if ((ctmp = bmap[cworkB[l]]) < cstart) {
1284307b7a18SHong Zhang                 *mat_j++ = ctmp;
1285307b7a18SHong Zhang                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1286307b7a18SHong Zhang                 mat_a   += bs2;
1287307b7a18SHong Zhang                 ilen_row++;
1288307b7a18SHong Zhang               } else break;
1289307b7a18SHong Zhang             }
1290307b7a18SHong Zhang             imark = l;
1291307b7a18SHong Zhang             for (l=0; l<nzA; l++) {
1292307b7a18SHong Zhang               *mat_j++ = cstart+cworkA[l];
1293307b7a18SHong Zhang               if (!ijonly) {
1294307b7a18SHong Zhang                 ierr   = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1295307b7a18SHong Zhang                 mat_a += bs2;
1296307b7a18SHong Zhang               }
1297307b7a18SHong Zhang               ilen_row++;
1298307b7a18SHong Zhang             }
1299307b7a18SHong Zhang             for (l=imark; l<nzB; l++) {
1300307b7a18SHong Zhang               *mat_j++ = bmap[cworkB[l]];
1301307b7a18SHong Zhang               if (!ijonly) {
1302307b7a18SHong Zhang                 ierr   = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1303307b7a18SHong Zhang                 mat_a += bs2;
1304307b7a18SHong Zhang               }
1305307b7a18SHong Zhang               ilen_row++;
1306307b7a18SHong Zhang             }
1307307b7a18SHong Zhang           }
1308d5b485abSSatish Balay           imat_ilen[row] = ilen_row;
1309d5b485abSSatish Balay         }
1310d5b485abSSatish Balay       }
1311d5b485abSSatish Balay     }
1312d5b485abSSatish Balay   }
1313d5b485abSSatish Balay 
1314d5b485abSSatish Balay   /*   Now assemble the off proc rows*/
1315d5b485abSSatish Balay   {
1316b24ad042SBarry Smith     PetscInt    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1317b24ad042SBarry Smith     PetscInt    *imat_j,*imat_i;
13180298fd71SBarry Smith     MatScalar   *imat_a = NULL,*rbuf4_i = NULL;
1319b24ad042SBarry Smith     PetscMPIInt ii;
1320d5b485abSSatish Balay 
1321d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
132226fbe8dcSKarl Rupp       if (ijonly) ii = tmp2;
132326fbe8dcSKarl Rupp       else {
1324b24ad042SBarry Smith         ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr);
13257a868f3eSHong Zhang       }
1326b24ad042SBarry Smith       idex    = pa[ii];
1327999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1328d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1329d5b485abSSatish Balay       ct1     = 2*jmax + 1;
1330d5b485abSSatish Balay       ct2     = 0;
1331b24ad042SBarry Smith       rbuf2_i = rbuf2[ii];
1332b24ad042SBarry Smith       rbuf3_i = rbuf3[ii];
13337a868f3eSHong Zhang       if (!ijonly) rbuf4_i = rbuf4[ii];
1334d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1335d5b485abSSatish Balay         is_no = sbuf1_i[2*j-1];
13368fa52d88SHong Zhang         if (!allcolumns[is_no]) cmap_i = cmap[is_no];
1337d5b485abSSatish Balay         rmap_i    = rmap[is_no];
1338df36731bSBarry Smith         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1339d5b485abSSatish Balay         imat_ilen = mat->ilen;
1340d5b485abSSatish Balay         imat_j    = mat->j;
1341d5b485abSSatish Balay         imat_i    = mat->i;
13427a868f3eSHong Zhang         if (!ijonly) imat_a = mat->a;
1343d5b485abSSatish Balay         max1 = sbuf1_i[2*j];
1344d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1345d5b485abSSatish Balay           row = sbuf1_i[ct1];
1346233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
13478fa52d88SHong Zhang           ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr);
1348233c2ff6SSatish Balay           row--;
134926fbe8dcSKarl Rupp           if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1350233c2ff6SSatish Balay #else
1351d5b485abSSatish Balay           row = rmap_i[row];
1352233c2ff6SSatish Balay #endif
1353d5b485abSSatish Balay           ilen  = imat_ilen[row];
1354df36731bSBarry Smith           mat_i = imat_i[row];
13557a868f3eSHong Zhang           if (!ijonly) mat_a = imat_a + mat_i*bs2;
1356d5b485abSSatish Balay           mat_j = imat_j + mat_i;
1357d5b485abSSatish Balay           max2  = rbuf2_i[ct1];
1358307b7a18SHong Zhang 
1359307b7a18SHong Zhang           if (!allcolumns[is_no]) {
1360d5b485abSSatish Balay             for (l=0; l<max2; l++,ct2++) {
1361233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
13628fa52d88SHong Zhang               ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1363233c2ff6SSatish Balay               if (tcol) {
1364233c2ff6SSatish Balay #else
1365d5b485abSSatish Balay               if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1366233c2ff6SSatish Balay #endif
1367df36731bSBarry Smith                 *mat_j++ = tcol - 1;
13687a868f3eSHong Zhang                 if (!ijonly) {
13693eda8832SBarry Smith                   ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1370549d3d68SSatish Balay                   mat_a += bs2;
13717a868f3eSHong Zhang                 }
1372d5b485abSSatish Balay                 ilen++;
1373d5b485abSSatish Balay               }
1374d5b485abSSatish Balay             }
1375307b7a18SHong Zhang           } else { /* allcolumns */
1376307b7a18SHong Zhang             for (l=0; l<max2; l++,ct2++) {
1377307b7a18SHong Zhang               *mat_j++ = rbuf3_i[ct2];
1378307b7a18SHong Zhang               if (!ijonly) {
1379307b7a18SHong Zhang                 ierr   = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1380307b7a18SHong Zhang                 mat_a += bs2;
1381307b7a18SHong Zhang               }
1382307b7a18SHong Zhang               ilen++;
1383307b7a18SHong Zhang             }
1384307b7a18SHong Zhang           }
1385d5b485abSSatish Balay           imat_ilen[row] = ilen;
1386d5b485abSSatish Balay         }
1387d5b485abSSatish Balay       }
1388d5b485abSSatish Balay     }
1389d5b485abSSatish Balay   }
13907a868f3eSHong Zhang   if (!ijonly) {
1391606d414cSSatish Balay     ierr = PetscFree(r_status4);CHKERRQ(ierr);
1392606d414cSSatish Balay     ierr = PetscFree(r_waits4);CHKERRQ(ierr);
13930c468ba9SBarry Smith     if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);}
1394606d414cSSatish Balay     ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1395606d414cSSatish Balay     ierr = PetscFree(s_status4);CHKERRQ(ierr);
13967a868f3eSHong Zhang   }
1397d5b485abSSatish Balay 
1398d5b485abSSatish Balay   /* Restore the indices */
1399d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
14004da72fa9SHong Zhang     if (!allrows[i]) {
1401d5b485abSSatish Balay       ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
14024da72fa9SHong Zhang     }
1403307b7a18SHong Zhang     if (!allcolumns[i]) {
1404d5b485abSSatish Balay       ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1405d5b485abSSatish Balay     }
1406307b7a18SHong Zhang   }
1407d5b485abSSatish Balay 
1408d5b485abSSatish Balay   /* Destroy allocated memory */
140923bdbc58SBarry Smith #if defined(PETSC_USE_CTABLE)
141023bdbc58SBarry Smith   ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr);
141123bdbc58SBarry Smith #else
141223bdbc58SBarry Smith   ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr);
141323bdbc58SBarry Smith #endif
141423bdbc58SBarry Smith   ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr);
1415606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
1416d5b485abSSatish Balay 
141723bdbc58SBarry Smith   ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr);
1418606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1419606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1420d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
1421606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1422d5b485abSSatish Balay   }
1423d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1424606d414cSSatish Balay     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1425d5b485abSSatish Balay   }
142623bdbc58SBarry Smith   ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr);
1427606d414cSSatish Balay   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1428606d414cSSatish Balay   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1429606d414cSSatish Balay   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
14307a868f3eSHong Zhang   if (!ijonly) {
14317a868f3eSHong Zhang     for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);}
14327a868f3eSHong Zhang     ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1433606d414cSSatish Balay     ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1434606d414cSSatish Balay     ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
14357a868f3eSHong Zhang   }
1436d5b485abSSatish Balay 
1437233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE)
1438ff0794f7SSatish Balay   for (i=0; i<ismax; i++) {
14398fa52d88SHong Zhang     ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr);
1440233c2ff6SSatish Balay   }
1441233c2ff6SSatish Balay #endif
14428fa52d88SHong Zhang   ierr = PetscFree(rmap);CHKERRQ(ierr);
14438fa52d88SHong Zhang 
14448fa52d88SHong Zhang   for (i=0; i<ismax; i++) {
14458fa52d88SHong Zhang     if (!allcolumns[i]) {
14468fa52d88SHong Zhang #if defined(PETSC_USE_CTABLE)
14478fa52d88SHong Zhang       ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr);
14488fa52d88SHong Zhang #else
14498fa52d88SHong Zhang       ierr = PetscFree(cmap[i]);CHKERRQ(ierr);
14508fa52d88SHong Zhang #endif
14518fa52d88SHong Zhang     }
14528fa52d88SHong Zhang   }
14538fa52d88SHong Zhang   ierr = PetscFree(cmap);CHKERRQ(ierr);
1454606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1455d5b485abSSatish Balay 
1456d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
145736f4e84dSSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
145836f4e84dSSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1459d5b485abSSatish Balay   }
14607a868f3eSHong Zhang 
14617a868f3eSHong Zhang   c->ijonly = PETSC_FALSE; /* set back to the default */
14623a40ed3dSBarry Smith   PetscFunctionReturn(0);
1463d5b485abSSatish Balay }
1464ca161407SBarry Smith 
1465