xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1d5b485abSSatish Balay /*
2d5b485abSSatish Balay    Routines to compute overlapping regions of a parallel MPI matrix
3d5b485abSSatish Balay   and to find submatrices that were shared across processors.
4d5b485abSSatish Balay */
570f55243SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h"
60a835dfdSSatish Balay #include "petscbt.h"
7d5b485abSSatish Balay 
8*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat,int,IS *);
9*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,int,char **,int*,int**);
10*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,int,int **,int**,int*);
11dfbe8321SBarry Smith EXTERN PetscErrorCode MatGetRow_MPIBAIJ(Mat,int,int*,int**,PetscScalar**);
12dfbe8321SBarry Smith EXTERN PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,int,int*,int**,PetscScalar**);
13d5b485abSSatish Balay 
144a2ae208SSatish Balay #undef __FUNCT__
154a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ"
16dfbe8321SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,int imax,IS is[],int ov)
17d5b485abSSatish Balay {
18d9489beaSHong Zhang   Mat_MPIBAIJ  *c = (Mat_MPIBAIJ*)C->data;
19*6849ba73SBarry Smith   PetscErrorCode ierr;
20*6849ba73SBarry Smith   int          i,N=C->N, bs=c->bs;
2136f4e84dSSatish Balay   IS           *is_new;
2236f4e84dSSatish Balay 
233a40ed3dSBarry Smith   PetscFunctionBegin;
2482502324SSatish Balay   ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr);
25df36731bSBarry Smith   /* Convert the indices into block format */
2627f478b1SHong Zhang   ierr = ISCompressIndicesGeneral(N,bs,imax,is,is_new);CHKERRQ(ierr);
2729bbc08cSBarry Smith   if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
28d5b485abSSatish Balay   for (i=0; i<ov; ++i) {
2936f4e84dSSatish Balay     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr);
30d5b485abSSatish Balay   }
31ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);}
3227f478b1SHong Zhang   ierr = ISExpandIndicesGeneral(N,bs,imax,is_new,is);CHKERRQ(ierr);
33ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);}
34606d414cSSatish Balay   ierr = PetscFree(is_new);CHKERRQ(ierr);
353a40ed3dSBarry Smith   PetscFunctionReturn(0);
36d5b485abSSatish Balay }
37d5b485abSSatish Balay 
38d5b485abSSatish Balay /*
39d5b485abSSatish Balay   Sample message format:
40d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
410e9b0e7eSHong Zhang   to index sets is[1], is[5]
42d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
43d5b485abSSatish Balay   -----------
44d5b485abSSatish Balay   mesg [1] = 1 => is[1]
45d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
46d5b485abSSatish Balay   -----------
47d5b485abSSatish Balay   mesg [5] = 5  => is[5]
48d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
49d5b485abSSatish Balay   -----------
50d5b485abSSatish Balay   mesg [7]
510e9b0e7eSHong Zhang   mesg [n]  data(is[1])
52d5b485abSSatish Balay   -----------
53d5b485abSSatish Balay   mesg[n+1]
54d5b485abSSatish Balay   mesg[m]  data(is[5])
55d5b485abSSatish Balay   -----------
56d5b485abSSatish Balay 
57d5b485abSSatish Balay   Notes:
58d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
59d5b485abSSatish Balay   nrqr - no of requests recieved (which have to be or which have been processed
60d5b485abSSatish Balay */
614a2ae208SSatish Balay #undef __FUNCT__
624a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once"
63*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,int imax,IS is[])
64d5b485abSSatish Balay {
65df36731bSBarry Smith   Mat_MPIBAIJ  *c = (Mat_MPIBAIJ*)C->data;
66d5b485abSSatish Balay   int         **idx,*n,*w1,*w2,*w3,*w4,*rtable,**data,len,*idx_i;
67*6849ba73SBarry Smith   PetscErrorCode ierr;
68*6849ba73SBarry Smith   int         size,rank,Mbs,i,j,k,**rbuf,row,proc,nrqs,msz,**outdat,**ptr;
69c7dd2924SSatish Balay   int         *ctr,*pa,*tmp,nrqr,*isz,*isz1,**xdata,**rbuf2;
70c7dd2924SSatish Balay   int         *onodes1,*olengths1,tag1,tag2,*onodes2,*olengths2;
71f1af5d2fSBarry Smith   PetscBT     *table;
72d5b485abSSatish Balay   MPI_Comm    comm;
73d5b485abSSatish Balay   MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2;
74d5b485abSSatish Balay   MPI_Status  *s_status,*recv_status;
75d5b485abSSatish Balay 
763a40ed3dSBarry Smith   PetscFunctionBegin;
77d5b485abSSatish Balay   comm   = C->comm;
78d5b485abSSatish Balay   size   = c->size;
79d5b485abSSatish Balay   rank   = c->rank;
80df36731bSBarry Smith   Mbs    = c->Mbs;
81d5b485abSSatish Balay 
82c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
83c7dd2924SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
84c7dd2924SSatish Balay 
85df36731bSBarry Smith   len    = (imax+1)*sizeof(int*)+ (imax + Mbs)*sizeof(int);
8682502324SSatish Balay   ierr   = PetscMalloc(len,&idx);CHKERRQ(ierr);
87d5b485abSSatish Balay   n      = (int*)(idx + imax);
88d5b485abSSatish Balay   rtable = n + imax;
89d5b485abSSatish Balay 
90d5b485abSSatish Balay   for (i=0; i<imax; i++) {
91d5b485abSSatish Balay     ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr);
92b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
93d5b485abSSatish Balay   }
94d5b485abSSatish Balay 
95d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
96d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
97d5b485abSSatish Balay     len = c->rowners[i+1];
98d5b485abSSatish Balay     for (; j<len; j++) {
99d5b485abSSatish Balay       rtable[j] = i;
100d5b485abSSatish Balay     }
101d5b485abSSatish Balay   }
102d5b485abSSatish Balay 
103d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
104d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
10582502324SSatish Balay   ierr = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr);/*  mesg size */
106d5b485abSSatish Balay   w2   = w1 + size;       /* if w2[i] marked, then a message to proc i*/
107d5b485abSSatish Balay   w3   = w2 + size;       /* no of IS that needs to be sent to proc i */
108d5b485abSSatish Balay   w4   = w3 + size;       /* temp work space used in determining w1, w2, w3 */
109549d3d68SSatish Balay   ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/
110d5b485abSSatish Balay   for (i=0; i<imax; i++) {
111549d3d68SSatish Balay     ierr  = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/
112d5b485abSSatish Balay     idx_i = idx[i];
113d5b485abSSatish Balay     len   = n[i];
114d5b485abSSatish Balay     for (j=0; j<len; j++) {
115d5b485abSSatish Balay       row  = idx_i[j];
1166b41c64dSBarry Smith       if (row < 0) {
11729bbc08cSBarry Smith         SETERRQ(1,"Index set cannot have negative entries");
1186b41c64dSBarry Smith       }
119d5b485abSSatish Balay       proc = rtable[row];
120d5b485abSSatish Balay       w4[proc]++;
121d5b485abSSatish Balay     }
122d5b485abSSatish Balay     for (j=0; j<size; j++){
123d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
124d5b485abSSatish Balay     }
125d5b485abSSatish Balay   }
126d5b485abSSatish Balay 
127d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
128d5b485abSSatish Balay   msz      = 0;              /* total mesg length (for all proc */
1290e9b0e7eSHong Zhang   w1[rank] = 0;              /* no mesg sent to itself */
130d5b485abSSatish Balay   w3[rank] = 0;
131d5b485abSSatish Balay   for (i=0; i<size; i++) {
132d5b485abSSatish Balay     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
133d5b485abSSatish Balay   }
134d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
135b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr);
136d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
137d5b485abSSatish Balay     if (w1[i]) {pa[j] = i; j++;}
138d5b485abSSatish Balay   }
139d5b485abSSatish Balay 
140d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
141d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
142d5b485abSSatish Balay     j      = pa[i];
143d5b485abSSatish Balay     w1[j] += w2[j] + 2*w3[j];
144d5b485abSSatish Balay     msz   += w1[j];
145d5b485abSSatish Balay   }
146d5b485abSSatish Balay 
147c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
148a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
149c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
150d5b485abSSatish Balay 
151c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
152c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr);
153d5b485abSSatish Balay 
154d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
155d5b485abSSatish Balay   len    = 2*size*sizeof(int*) + (size+msz)*sizeof(int);
15682502324SSatish Balay   ierr = PetscMalloc(len,&outdat);CHKERRQ(ierr);
157d5b485abSSatish Balay   ptr    = outdat + size;     /* Pointers to the data in outgoing buffers */
158549d3d68SSatish Balay   ierr   = PetscMemzero(outdat,2*size*sizeof(int*));CHKERRQ(ierr);
159d5b485abSSatish Balay   tmp    = (int*)(outdat + 2*size);
160d5b485abSSatish Balay   ctr    = tmp + msz;
161d5b485abSSatish Balay 
162d5b485abSSatish Balay   {
163d5b485abSSatish Balay     int *iptr = tmp,ict  = 0;
164d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
165d5b485abSSatish Balay       j         = pa[i];
166d5b485abSSatish Balay       iptr     +=  ict;
167d5b485abSSatish Balay       outdat[j] = iptr;
168d5b485abSSatish Balay       ict       = w1[j];
169d5b485abSSatish Balay     }
170d5b485abSSatish Balay   }
171d5b485abSSatish Balay 
172d5b485abSSatish Balay   /* Form the outgoing messages */
173d5b485abSSatish Balay   /*plug in the headers*/
174d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
175d5b485abSSatish Balay     j            = pa[i];
176d5b485abSSatish Balay     outdat[j][0] = 0;
177549d3d68SSatish Balay     ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr);
178d5b485abSSatish Balay     ptr[j]       = outdat[j] + 2*w3[j] + 1;
179d5b485abSSatish Balay   }
180d5b485abSSatish Balay 
181d5b485abSSatish Balay   /* Memory for doing local proc's work*/
182d5b485abSSatish Balay   {
183d5b485abSSatish Balay     int  *d_p;
184d5b485abSSatish Balay     char *t_p;
185d5b485abSSatish Balay 
186f1af5d2fSBarry Smith     len  = (imax)*(sizeof(PetscBT) + sizeof(int*)+ sizeof(int)) +
187b6410449SSatish Balay       (Mbs)*imax*sizeof(int)  + (Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char) + 1;
18882502324SSatish Balay     ierr = PetscMalloc(len,&table);CHKERRQ(ierr);
189549d3d68SSatish Balay     ierr = PetscMemzero(table,len);CHKERRQ(ierr);
190d5b485abSSatish Balay     data = (int **)(table + imax);
191bbd702dbSSatish Balay     isz  = (int  *)(data  + imax);
192bbd702dbSSatish Balay     d_p  = (int  *)(isz   + imax);
193bbd702dbSSatish Balay     t_p  = (char *)(d_p   + Mbs*imax);
194bbd702dbSSatish Balay     for (i=0; i<imax; i++) {
195b6410449SSatish Balay       table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
196bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
197d5b485abSSatish Balay     }
198d5b485abSSatish Balay   }
199d5b485abSSatish Balay 
200d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
201d5b485abSSatish Balay   {
202d5b485abSSatish Balay     int     n_i,*data_i,isz_i,*outdat_j,ctr_j;
203f1af5d2fSBarry Smith     PetscBT table_i;
204d5b485abSSatish Balay 
205d5b485abSSatish Balay     for (i=0; i<imax; i++) {
206549d3d68SSatish Balay       ierr    = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr);
207d5b485abSSatish Balay       n_i     = n[i];
208d5b485abSSatish Balay       table_i = table[i];
209d5b485abSSatish Balay       idx_i   = idx[i];
210d5b485abSSatish Balay       data_i  = data[i];
211d5b485abSSatish Balay       isz_i   = isz[i];
212d5b485abSSatish Balay       for (j=0;  j<n_i; j++) {  /* parse the indices of each IS */
213d5b485abSSatish Balay         row  = idx_i[j];
214d5b485abSSatish Balay         proc = rtable[row];
215d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
216d5b485abSSatish Balay           ctr[proc]++;
217d5b485abSSatish Balay           *ptr[proc] = row;
218d5b485abSSatish Balay           ptr[proc]++;
219d5b485abSSatish Balay         }
220d5b485abSSatish Balay         else { /* Update the local table */
221f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
222d5b485abSSatish Balay         }
223d5b485abSSatish Balay       }
224d5b485abSSatish Balay       /* Update the headers for the current IS */
225d5b485abSSatish Balay       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
226d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
227d5b485abSSatish Balay           outdat_j        = outdat[j];
228d5b485abSSatish Balay           k               = ++outdat_j[0];
229d5b485abSSatish Balay           outdat_j[2*k]   = ctr_j;
230d5b485abSSatish Balay           outdat_j[2*k-1] = i;
231d5b485abSSatish Balay         }
232d5b485abSSatish Balay       }
233d5b485abSSatish Balay       isz[i] = isz_i;
234d5b485abSSatish Balay     }
235d5b485abSSatish Balay   }
236d5b485abSSatish Balay 
237d5b485abSSatish Balay   /*  Now  post the sends */
238b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
239d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
240d5b485abSSatish Balay     j    = pa[i];
241c7dd2924SSatish Balay     ierr = MPI_Isend(outdat[j],w1[j],MPI_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr);
242d5b485abSSatish Balay   }
243d5b485abSSatish Balay 
244d5b485abSSatish Balay   /* No longer need the original indices*/
245d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
246d5b485abSSatish Balay     ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr);
247d5b485abSSatish Balay   }
248606d414cSSatish Balay   ierr = PetscFree(idx);CHKERRQ(ierr);
249d5b485abSSatish Balay 
250d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
251d5b485abSSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
252d5b485abSSatish Balay   }
253d5b485abSSatish Balay 
254d5b485abSSatish Balay   /* Do Local work*/
255df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
256d5b485abSSatish Balay 
257d5b485abSSatish Balay   /* Receive messages*/
258b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr);
259c7dd2924SSatish Balay   ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);
260d5b485abSSatish Balay 
261b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr);
262ca161407SBarry Smith   ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);
263d5b485abSSatish Balay 
264d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
265606d414cSSatish Balay   ierr = PetscFree(outdat);CHKERRQ(ierr);
266606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
267d5b485abSSatish Balay 
26882502324SSatish Balay   ierr = PetscMalloc((nrqr+1)*sizeof(int*),&xdata);CHKERRQ(ierr);
269b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(int),&isz1);CHKERRQ(ierr);
270df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
271606d414cSSatish Balay   ierr = PetscFree(rbuf);CHKERRQ(ierr);
272d5b485abSSatish Balay 
273d5b485abSSatish Balay   /* Send the data back*/
274d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
275d5b485abSSatish Balay   {
276c7dd2924SSatish Balay     int *rw1;
277d5b485abSSatish Balay 
278c7dd2924SSatish Balay     ierr = PetscMalloc(size*sizeof(int),&rw1);CHKERRQ(ierr);
279c7dd2924SSatish Balay     ierr = PetscMemzero(rw1,size*sizeof(int));CHKERRQ(ierr);
280c7dd2924SSatish Balay 
281d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
282d5b485abSSatish Balay       proc      = recv_status[i].MPI_SOURCE;
283c7dd2924SSatish Balay       if (proc != onodes1[i]) SETERRQ(1,"MPI_SOURCE mismatch");
284d5b485abSSatish Balay       rw1[proc] = isz1[i];
285d5b485abSSatish Balay     }
286d5b485abSSatish Balay 
287c7dd2924SSatish Balay     ierr = PetscFree(onodes1);CHKERRQ(ierr);
288c7dd2924SSatish Balay     ierr = PetscFree(olengths1);CHKERRQ(ierr);
289d5b485abSSatish Balay 
290c7dd2924SSatish Balay     /* Determine the number of messages to expect, their lengths, from from-ids */
291c7dd2924SSatish Balay     ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr);
292c7dd2924SSatish Balay     PetscFree(rw1);
293c7dd2924SSatish Balay   }
294c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
295c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr);
296d5b485abSSatish Balay 
297d5b485abSSatish Balay   /*  Now  post the sends */
298b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
299d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
300d5b485abSSatish Balay     j    = recv_status[i].MPI_SOURCE;
301c7dd2924SSatish Balay     ierr = MPI_Isend(xdata[i],isz1[i],MPI_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr);
302d5b485abSSatish Balay   }
303d5b485abSSatish Balay 
304d5b485abSSatish Balay   /* receive work done on other processors*/
305d5b485abSSatish Balay   {
306999d9058SBarry Smith     int         idex,is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
307f1af5d2fSBarry Smith     PetscBT     table_i;
308d5b485abSSatish Balay     MPI_Status  *status2;
309d5b485abSSatish Balay 
310169449a0SSatish Balay     ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr);
311d5b485abSSatish Balay 
312d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
313999d9058SBarry Smith       ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr);
314d5b485abSSatish Balay       /* Process the message*/
315999d9058SBarry Smith       rbuf2_i = rbuf2[idex];
316d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
317999d9058SBarry Smith       jmax    = rbuf2[idex][0];
318d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
319d5b485abSSatish Balay         max     = rbuf2_i[2*j];
320d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
321d5b485abSSatish Balay         isz_i   = isz[is_no];
322d5b485abSSatish Balay         data_i  = data[is_no];
323d5b485abSSatish Balay         table_i = table[is_no];
324d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
325d5b485abSSatish Balay           row = rbuf2_i[ct1];
326f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
327d5b485abSSatish Balay         }
328d5b485abSSatish Balay         isz[is_no] = isz_i;
329d5b485abSSatish Balay       }
330d5b485abSSatish Balay     }
331ca161407SBarry Smith     ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);
332606d414cSSatish Balay     ierr = PetscFree(status2);CHKERRQ(ierr);
333d5b485abSSatish Balay   }
334d5b485abSSatish Balay 
335d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
336029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr);
337d5b485abSSatish Balay   }
338d5b485abSSatish Balay 
339c7dd2924SSatish Balay 
340c7dd2924SSatish Balay   ierr = PetscFree(onodes2);CHKERRQ(ierr);
341c7dd2924SSatish Balay   ierr = PetscFree(olengths2);CHKERRQ(ierr);
342c7dd2924SSatish Balay 
343606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
344606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
345606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
346606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
347606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
348606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
349606d414cSSatish Balay   ierr = PetscFree(table);CHKERRQ(ierr);
350606d414cSSatish Balay   ierr = PetscFree(s_status);CHKERRQ(ierr);
351606d414cSSatish Balay   ierr = PetscFree(recv_status);CHKERRQ(ierr);
352606d414cSSatish Balay   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
353606d414cSSatish Balay   ierr = PetscFree(xdata);CHKERRQ(ierr);
354606d414cSSatish Balay   ierr = PetscFree(isz1);CHKERRQ(ierr);
3553a40ed3dSBarry Smith   PetscFunctionReturn(0);
356d5b485abSSatish Balay }
357d5b485abSSatish Balay 
3584a2ae208SSatish Balay #undef __FUNCT__
3594a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local"
360d5b485abSSatish Balay /*
361df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
362d5b485abSSatish Balay        the work on the local processor.
363d5b485abSSatish Balay 
364d5b485abSSatish Balay      Inputs:
365df36731bSBarry Smith       C      - MAT_MPIBAIJ;
366d5b485abSSatish Balay       imax - total no of index sets processed at a time;
367df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
368d5b485abSSatish Balay 
369d5b485abSSatish Balay      Output:
3700e9b0e7eSHong Zhang       isz    - array containing the count of the solution elements corresponding
371d5b485abSSatish Balay                to each index set;
372d5b485abSSatish Balay       data   - pointer to the solutions
373d5b485abSSatish Balay */
374*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,int imax,PetscBT *table,int *isz,int **data)
375d5b485abSSatish Balay {
376df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
377d5b485abSSatish Balay   Mat         A = c->A,B = c->B;
378df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
379df36731bSBarry Smith   int         start,end,val,max,rstart,cstart,*ai,*aj;
380d5b485abSSatish Balay   int         *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
381f1af5d2fSBarry Smith   PetscBT     table_i;
382d5b485abSSatish Balay 
3833a40ed3dSBarry Smith   PetscFunctionBegin;
384d5b485abSSatish Balay   rstart = c->rstart;
385d5b485abSSatish Balay   cstart = c->cstart;
386d5b485abSSatish Balay   ai     = a->i;
387df36731bSBarry Smith   aj     = a->j;
388d5b485abSSatish Balay   bi     = b->i;
389df36731bSBarry Smith   bj     = b->j;
390d5b485abSSatish Balay   garray = c->garray;
391d5b485abSSatish Balay 
392d5b485abSSatish Balay 
393d5b485abSSatish Balay   for (i=0; i<imax; i++) {
394d5b485abSSatish Balay     data_i  = data[i];
395d5b485abSSatish Balay     table_i = table[i];
396d5b485abSSatish Balay     isz_i   = isz[i];
397d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
398d5b485abSSatish Balay       row   = data_i[j] - rstart;
399d5b485abSSatish Balay       start = ai[row];
400d5b485abSSatish Balay       end   = ai[row+1];
401d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
402df36731bSBarry Smith         val = aj[k] + cstart;
403f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
404d5b485abSSatish Balay       }
405d5b485abSSatish Balay       start = bi[row];
406d5b485abSSatish Balay       end   = bi[row+1];
407d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
408df36731bSBarry Smith         val = garray[bj[k]];
409f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
410d5b485abSSatish Balay       }
411d5b485abSSatish Balay     }
412d5b485abSSatish Balay     isz[i] = isz_i;
413d5b485abSSatish Balay   }
4143a40ed3dSBarry Smith   PetscFunctionReturn(0);
415d5b485abSSatish Balay }
4164a2ae208SSatish Balay #undef __FUNCT__
4174a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive"
418d5b485abSSatish Balay /*
419df36731bSBarry Smith       MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages,
420d5b485abSSatish Balay          and return the output
421d5b485abSSatish Balay 
422d5b485abSSatish Balay          Input:
423d5b485abSSatish Balay            C    - the matrix
424d5b485abSSatish Balay            nrqr - no of messages being processed.
425d5b485abSSatish Balay            rbuf - an array of pointers to the recieved requests
426d5b485abSSatish Balay 
427d5b485abSSatish Balay          Output:
428d5b485abSSatish Balay            xdata - array of messages to be sent back
429d5b485abSSatish Balay            isz1  - size of each message
430d5b485abSSatish Balay 
431d5b485abSSatish Balay   For better efficiency perhaps we should malloc seperately each xdata[i],
432d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
4330e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of
434d5b485abSSatish Balay memory is used.
435d5b485abSSatish Balay 
436d5b485abSSatish Balay */
437*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,int nrqr,int **rbuf,int **xdata,int * isz1)
438d5b485abSSatish Balay {
439df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
440d5b485abSSatish Balay   Mat         A = c->A,B = c->B;
441df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
442*6849ba73SBarry Smith   PetscErrorCode ierr;
443df36731bSBarry Smith   int         rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
444d5b485abSSatish Balay   int         row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
445df36731bSBarry Smith   int         val,max1,max2,rank,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
446*6849ba73SBarry Smith   int         *rbuf_i,kmax,rbuf_0;
447f1af5d2fSBarry Smith   PetscBT     xtable;
448d5b485abSSatish Balay 
4493a40ed3dSBarry Smith   PetscFunctionBegin;
450d5b485abSSatish Balay   rank   = c->rank;
451df36731bSBarry Smith   Mbs    = c->Mbs;
452d5b485abSSatish Balay   rstart = c->rstart;
453d5b485abSSatish Balay   cstart = c->cstart;
454d5b485abSSatish Balay   ai     = a->i;
455df36731bSBarry Smith   aj     = a->j;
456d5b485abSSatish Balay   bi     = b->i;
457df36731bSBarry Smith   bj     = b->j;
458d5b485abSSatish Balay   garray = c->garray;
459d5b485abSSatish Balay 
460d5b485abSSatish Balay 
461d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
462d5b485abSSatish Balay     rbuf_i  =  rbuf[i];
463d5b485abSSatish Balay     rbuf_0  =  rbuf_i[0];
464d5b485abSSatish Balay     ct     += rbuf_0;
465d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; }
466d5b485abSSatish Balay   }
467d5b485abSSatish Balay 
468701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
469701b8814SBarry Smith   else        max1 = 1;
470d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
47182502324SSatish Balay   ierr         = PetscMalloc(mem_estimate*sizeof(int),&xdata[0]);CHKERRQ(ierr);
472d5b485abSSatish Balay   ++no_malloc;
4736831982aSBarry Smith   ierr         = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr);
474549d3d68SSatish Balay   ierr         = PetscMemzero(isz1,nrqr*sizeof(int));CHKERRQ(ierr);
475d5b485abSSatish Balay 
476d5b485abSSatish Balay   ct3 = 0;
477d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
478d5b485abSSatish Balay     rbuf_i =  rbuf[i];
479d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
480d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
481d5b485abSSatish Balay     ct2    =  ct1;
482d5b485abSSatish Balay     ct3    += ct1;
483d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
4846831982aSBarry Smith       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
485d5b485abSSatish Balay       oct2 = ct2;
486d5b485abSSatish Balay       kmax = rbuf_i[2*j];
487d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
488d5b485abSSatish Balay         row = rbuf_i[ct1];
489f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
490d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
491d5b485abSSatish Balay             new_estimate = (int)(1.5*mem_estimate)+1;
49282502324SSatish Balay             ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr);
493549d3d68SSatish Balay             ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr);
494606d414cSSatish Balay             ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
495d5b485abSSatish Balay             xdata[0]     = tmp;
496d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
497d5b485abSSatish Balay             for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
498d5b485abSSatish Balay           }
499d5b485abSSatish Balay           xdata[i][ct2++] = row;
500d5b485abSSatish Balay           ct3++;
501d5b485abSSatish Balay         }
502d5b485abSSatish Balay       }
503d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
504d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
505d5b485abSSatish Balay         start = ai[row];
506d5b485abSSatish Balay         end   = ai[row+1];
507d5b485abSSatish Balay         for (l=start; l<end; l++) {
508df36731bSBarry Smith           val = aj[l] + cstart;
509f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
510d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
511d5b485abSSatish Balay               new_estimate = (int)(1.5*mem_estimate)+1;
51282502324SSatish Balay               ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr);
513549d3d68SSatish Balay               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr);
514606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
515d5b485abSSatish Balay               xdata[0]     = tmp;
516d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
517d5b485abSSatish Balay               for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
518d5b485abSSatish Balay             }
519d5b485abSSatish Balay             xdata[i][ct2++] = val;
520d5b485abSSatish Balay             ct3++;
521d5b485abSSatish Balay           }
522d5b485abSSatish Balay         }
523d5b485abSSatish Balay         start = bi[row];
524d5b485abSSatish Balay         end   = bi[row+1];
525d5b485abSSatish Balay         for (l=start; l<end; l++) {
526df36731bSBarry Smith           val = garray[bj[l]];
527f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
528d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
529d5b485abSSatish Balay               new_estimate = (int)(1.5*mem_estimate)+1;
53082502324SSatish Balay               ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr);
531549d3d68SSatish Balay               ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr);
532606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
533d5b485abSSatish Balay               xdata[0]     = tmp;
534d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
535d5b485abSSatish Balay               for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
536d5b485abSSatish Balay             }
537d5b485abSSatish Balay             xdata[i][ct2++] = val;
538d5b485abSSatish Balay             ct3++;
539d5b485abSSatish Balay           }
540d5b485abSSatish Balay         }
541d5b485abSSatish Balay       }
542d5b485abSSatish Balay       /* Update the header*/
543d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
544d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
545d5b485abSSatish Balay     }
546d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
547d5b485abSSatish Balay     xdata[i+1]  = xdata[i] + ct2;
548d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
549d5b485abSSatish Balay   }
5506831982aSBarry Smith   ierr = PetscBTDestroy(xtable);CHKERRQ(ierr);
551b0a32e0cSBarry Smith   PetscLogInfo(0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %d bytes, required %d, no of mallocs = %d\n",rank,mem_estimate,ct3,no_malloc);
5523a40ed3dSBarry Smith   PetscFunctionReturn(0);
553d5b485abSSatish Balay }
554d5b485abSSatish Balay 
555*6849ba73SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat,int,const IS[],const IS[],MatReuse,Mat *);
556a2feddc5SSatish Balay 
5574a2ae208SSatish Balay #undef __FUNCT__
5584a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ"
559dfbe8321SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
560d5b485abSSatish Balay {
56136f4e84dSSatish Balay   IS          *isrow_new,*iscol_new;
562cf4f527aSSatish Balay   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
563*6849ba73SBarry Smith   PetscErrorCode ierr;
564*6849ba73SBarry Smith   int         nmax,nstages_local,nstages,i,pos,max_no,N=C->N,bs=c->bs;
565a2feddc5SSatish Balay 
5663a40ed3dSBarry Smith   PetscFunctionBegin;
5673a40ed3dSBarry Smith   /* The compression and expansion should be avoided. Does'nt point
568a2feddc5SSatish Balay      out errors might change the indices hence buggey */
569a2feddc5SSatish Balay 
570b0a32e0cSBarry Smith   ierr = PetscMalloc(2*(ismax+1)*sizeof(IS),&isrow_new);CHKERRQ(ierr);
57136f4e84dSSatish Balay   iscol_new = isrow_new + ismax;
57227f478b1SHong Zhang   ierr = ISCompressIndicesSorted(N,bs,ismax,isrow,isrow_new);CHKERRQ(ierr);
57327f478b1SHong Zhang   ierr = ISCompressIndicesSorted(N,bs,ismax,iscol,iscol_new);CHKERRQ(ierr);
574cf4f527aSSatish Balay 
575cf4f527aSSatish Balay   /* Allocate memory to hold all the submatrices */
576cf4f527aSSatish Balay   if (scall != MAT_REUSE_MATRIX) {
57782502324SSatish Balay     ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr);
578cf4f527aSSatish Balay   }
579cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
580cf4f527aSSatish Balay   nmax          = 20*1000000 / (c->Nbs * sizeof(int));
581329f5518SBarry Smith   if (!nmax) nmax = 1;
582cf4f527aSSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
583cf4f527aSSatish Balay 
584653e4784SBarry Smith   /* Make sure every processor loops through the nstages */
585ca161407SBarry Smith   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr);
586cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
587cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
588cf4f527aSSatish Balay     else if (pos == ismax) max_no = 0;
589cf4f527aSSatish Balay     else                   max_no = ismax-pos;
590cf4f527aSSatish Balay     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr);
591cf4f527aSSatish Balay     pos += max_no;
592cf4f527aSSatish Balay   }
59336f4e84dSSatish Balay 
59436f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
595ca161407SBarry Smith     ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr);
596ca161407SBarry Smith     ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr);
59736f4e84dSSatish Balay   }
598606d414cSSatish Balay   ierr = PetscFree(isrow_new);CHKERRQ(ierr);
5993a40ed3dSBarry Smith   PetscFunctionReturn(0);
600a2feddc5SSatish Balay }
601a2feddc5SSatish Balay 
602233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
6034a2ae208SSatish Balay #undef __FUNCT__
6044a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc"
605dfbe8321SBarry Smith PetscErrorCode PetscGetProc(const int gid, const int size, const int proc_gnode[], int *proc)
606233c2ff6SSatish Balay {
60723ce1328SBarry Smith   int nGlobalNd = proc_gnode[size];
60823ce1328SBarry Smith   int fproc = (int) ((float)gid * (float)size / (float)nGlobalNd + 0.5);
609233c2ff6SSatish Balay 
610233c2ff6SSatish Balay   PetscFunctionBegin;
61129bbc08cSBarry Smith   /* if(fproc < 0) SETERRQ(1,"fproc < 0");*/
61223ce1328SBarry Smith   if (fproc > size) fproc = size;
613233c2ff6SSatish Balay   while (gid < proc_gnode[fproc] || gid >= proc_gnode[fproc+1]) {
614233c2ff6SSatish Balay     if (gid < proc_gnode[fproc]) fproc--;
615233c2ff6SSatish Balay     else                         fproc++;
616233c2ff6SSatish Balay   }
61723ce1328SBarry Smith   /* if(fproc<0 || fproc>=size) { SETERRQ(1,"fproc < 0 || fproc >= size"); }*/
618233c2ff6SSatish Balay   *proc = fproc;
619233c2ff6SSatish Balay   PetscFunctionReturn(0);
620233c2ff6SSatish Balay }
621233c2ff6SSatish Balay #endif
622233c2ff6SSatish Balay 
623a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
6244a2ae208SSatish Balay #undef __FUNCT__
6254a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local"
626*6849ba73SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
627a2feddc5SSatish Balay {
628df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
629cf4f527aSSatish Balay   Mat         A = c->A;
630df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
631233c2ff6SSatish Balay   int         **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,start,end,size;
632*6849ba73SBarry Smith   PetscErrorCode ierr;
633*6849ba73SBarry Smith   int         **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row,proc;
634999d9058SBarry Smith   int         nrqs,msz,**ptr,idex,*req_size,*ctr,*pa,*tmp,tcol,nrqr;
635233c2ff6SSatish Balay   int         **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
636233c2ff6SSatish Balay   int         **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i;
637233c2ff6SSatish Balay   int         len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
638233c2ff6SSatish Balay   int         bs=c->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
63936f4e84dSSatish Balay   int         cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
64034e6ae18SSatish Balay   int         *bmap = c->garray,ctmp,rstart=c->rstart,tag0,tag1,tag2,tag3;
641d5b485abSSatish Balay   MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
642d5b485abSSatish Balay   MPI_Request *r_waits4,*s_waits3,*s_waits4;
643d5b485abSSatish Balay   MPI_Status  *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
644d5b485abSSatish Balay   MPI_Status  *r_status3,*r_status4,*s_status4;
645d5b485abSSatish Balay   MPI_Comm    comm;
6463eda8832SBarry Smith   MatScalar   **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB;
6473eda8832SBarry Smith   MatScalar   *a_a=a->a,*b_a=b->a;
6480f5bd95cSBarry Smith   PetscTruth  flag;
649c7dd2924SSatish Balay   int         *onodes1,*olengths1;
650c7dd2924SSatish Balay 
65180d608c0SSatish Balay #if defined (PETSC_USE_CTABLE)
652233c2ff6SSatish Balay   int tt;
653233c2ff6SSatish Balay   PetscTable  *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1;
654233c2ff6SSatish Balay #else
655233c2ff6SSatish Balay   int         **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
656233c2ff6SSatish Balay #endif
657d5b485abSSatish Balay 
6583a40ed3dSBarry Smith   PetscFunctionBegin;
659d5b485abSSatish Balay   comm   = C->comm;
66034e6ae18SSatish Balay   tag0   = C->tag;
661d5b485abSSatish Balay   size   = c->size;
662d5b485abSSatish Balay   rank   = c->rank;
663d5b485abSSatish Balay 
66434e6ae18SSatish Balay   /* Get some new tags to keep the communication clean */
66534e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
66634e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
66734e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
66834e6ae18SSatish Balay 
669d5b485abSSatish Balay   /* Check if the col indices are sorted */
670d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
671888f2ed8SSatish Balay     ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr);
67229bbc08cSBarry Smith     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
673d5b485abSSatish Balay   }
674d5b485abSSatish Balay 
675233c2ff6SSatish Balay   len    = (2*ismax+1)*(sizeof(int*)+ sizeof(int));
676233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE)
677233c2ff6SSatish Balay   len    += (Mbs+1)*sizeof(int);
678233c2ff6SSatish Balay #endif
67982502324SSatish Balay   ierr = PetscMalloc(len,&irow);CHKERRQ(ierr);
680d5b485abSSatish Balay   icol = irow + ismax;
681d5b485abSSatish Balay   nrow = (int*)(icol + ismax);
682d5b485abSSatish Balay   ncol = nrow + ismax;
683233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE)
684d5b485abSSatish Balay   rtable = ncol + ismax;
685d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
686d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
687d5b485abSSatish Balay     jmax = c->rowners[i+1];
688d5b485abSSatish Balay     for (; j<jmax; j++) {
689d5b485abSSatish Balay       rtable[j] = i;
690d5b485abSSatish Balay     }
691d5b485abSSatish Balay   }
692233c2ff6SSatish Balay #endif
693233c2ff6SSatish Balay 
694233c2ff6SSatish Balay   for (i=0; i<ismax; i++) {
695233c2ff6SSatish Balay     ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
696233c2ff6SSatish Balay     ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
697233c2ff6SSatish Balay     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
698233c2ff6SSatish Balay     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
699233c2ff6SSatish Balay   }
700d5b485abSSatish Balay 
701d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg,and buffer space
702d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
70382502324SSatish Balay   ierr = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr); /* mesg size */
704d5b485abSSatish Balay   w2   = w1 + size;      /* if w2[i] marked, then a message to proc i*/
705d5b485abSSatish Balay   w3   = w2 + size;      /* no of IS that needs to be sent to proc i */
706d5b485abSSatish Balay   w4   = w3 + size;      /* temp work space used in determining w1, w2, w3 */
707549d3d68SSatish Balay   ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/
708d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
709549d3d68SSatish Balay     ierr   = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/
710d5b485abSSatish Balay     jmax   = nrow[i];
711d5b485abSSatish Balay     irow_i = irow[i];
712d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
713d5b485abSSatish Balay       row  = irow_i[j];
714233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
715233c2ff6SSatish Balay       ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr);
716233c2ff6SSatish Balay #else
717d5b485abSSatish Balay       proc = rtable[row];
718233c2ff6SSatish Balay #endif
719d5b485abSSatish Balay       w4[proc]++;
720d5b485abSSatish Balay     }
721d5b485abSSatish Balay     for (j=0; j<size; j++) {
722d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
723d5b485abSSatish Balay     }
724d5b485abSSatish Balay   }
725d5b485abSSatish Balay 
726d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
727e757655aSSatish Balay   msz      = 0;              /* total mesg length for all proc */
728d5b485abSSatish Balay   w1[rank] = 0;              /* no mesg sent to intself */
729d5b485abSSatish Balay   w3[rank] = 0;
730d5b485abSSatish Balay   for (i=0; i<size; i++) {
731d5b485abSSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
732d5b485abSSatish Balay   }
733b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr); /*(proc -array)*/
734d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
735d5b485abSSatish Balay     if (w1[i]) { pa[j] = i; j++; }
736d5b485abSSatish Balay   }
737d5b485abSSatish Balay 
738d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
739d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
740d5b485abSSatish Balay     j     = pa[i];
741d5b485abSSatish Balay     w1[j] += w2[j] + 2* w3[j];
742d5b485abSSatish Balay     msz   += w1[j];
743d5b485abSSatish Balay   }
744d5b485abSSatish Balay 
745c7dd2924SSatish Balay   /* Determine the number of messages to expect, their lengths, from from-ids */
746a96d083dSSatish Balay   ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr);
747c7dd2924SSatish Balay   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr);
748d5b485abSSatish Balay 
749c7dd2924SSatish Balay   /* Now post the Irecvs corresponding to these messages */
750c7dd2924SSatish Balay   ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr);
751c7dd2924SSatish Balay 
752c7dd2924SSatish Balay   ierr = PetscFree(onodes1);CHKERRQ(ierr);
753c7dd2924SSatish Balay   ierr = PetscFree(olengths1);CHKERRQ(ierr);
754d5b485abSSatish Balay 
755d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
756d5b485abSSatish Balay   len  = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int);
75782502324SSatish Balay   ierr = PetscMalloc(len,&sbuf1);CHKERRQ(ierr);
758d5b485abSSatish Balay   ptr  = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
759549d3d68SSatish Balay   ierr = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr);
760d5b485abSSatish Balay   /* allocate memory for outgoing data + buf to receive the first reply */
761d5b485abSSatish Balay   tmp  = (int*)(ptr + size);
762d5b485abSSatish Balay   ctr  = tmp + 2*msz;
763d5b485abSSatish Balay 
764d5b485abSSatish Balay   {
765d5b485abSSatish Balay     int *iptr = tmp,ict = 0;
766d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
767d5b485abSSatish Balay       j         = pa[i];
768d5b485abSSatish Balay       iptr     += ict;
769d5b485abSSatish Balay       sbuf1[j]  = iptr;
770d5b485abSSatish Balay       ict       = w1[j];
771d5b485abSSatish Balay     }
772d5b485abSSatish Balay   }
773d5b485abSSatish Balay 
774d5b485abSSatish Balay   /* Form the outgoing messages */
775d5b485abSSatish Balay   /* Initialise the header space */
776d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
777d5b485abSSatish Balay     j           = pa[i];
778d5b485abSSatish Balay     sbuf1[j][0] = 0;
779549d3d68SSatish Balay     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr);
780d5b485abSSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
781d5b485abSSatish Balay   }
782d5b485abSSatish Balay 
783d5b485abSSatish Balay   /* Parse the isrow and copy data into outbuf */
784d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
785549d3d68SSatish Balay     ierr   = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr);
786d5b485abSSatish Balay     irow_i = irow[i];
787d5b485abSSatish Balay     jmax   = nrow[i];
788d5b485abSSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
789d5b485abSSatish Balay       row  = irow_i[j];
790233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
791233c2ff6SSatish Balay       ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr);
792233c2ff6SSatish Balay #else
793d5b485abSSatish Balay       proc = rtable[row];
794233c2ff6SSatish Balay #endif
795d5b485abSSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
796d5b485abSSatish Balay         ctr[proc]++;
797d5b485abSSatish Balay         *ptr[proc] = row;
798d5b485abSSatish Balay         ptr[proc]++;
799d5b485abSSatish Balay       }
800d5b485abSSatish Balay     }
801d5b485abSSatish Balay     /* Update the headers for the current IS */
802d5b485abSSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
80306ef35abSSatish Balay       if ((ctr_j = ctr[j])) {
804d5b485abSSatish Balay         sbuf1_j        = sbuf1[j];
805d5b485abSSatish Balay         k              = ++sbuf1_j[0];
806d5b485abSSatish Balay         sbuf1_j[2*k]   = ctr_j;
807d5b485abSSatish Balay         sbuf1_j[2*k-1] = i;
808d5b485abSSatish Balay       }
809d5b485abSSatish Balay     }
810d5b485abSSatish Balay   }
811d5b485abSSatish Balay 
812d5b485abSSatish Balay   /*  Now  post the sends */
813b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
814d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
815d5b485abSSatish Balay     j = pa[i];
816ca161407SBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
817d5b485abSSatish Balay   }
818d5b485abSSatish Balay 
819d5b485abSSatish Balay   /* Post Recieves to capture the buffer size */
820b0a32e0cSBarry Smith   ierr     = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr);
82182502324SSatish Balay   ierr     = PetscMalloc((nrqs+1)*sizeof(int*),&rbuf2);CHKERRQ(ierr);
822d5b485abSSatish Balay   rbuf2[0] = tmp + msz;
823d5b485abSSatish Balay   for (i=1; i<nrqs; ++i) {
824d5b485abSSatish Balay     j        = pa[i];
825d5b485abSSatish Balay     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
826d5b485abSSatish Balay   }
827d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
828d5b485abSSatish Balay     j    = pa[i];
829ca161407SBarry Smith     ierr = MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
830d5b485abSSatish Balay   }
831d5b485abSSatish Balay 
832d5b485abSSatish Balay   /* Send to other procs the buf size they should allocate */
833d5b485abSSatish Balay 
834d5b485abSSatish Balay   /* Receive messages*/
835b0a32e0cSBarry Smith   ierr       = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
836b0a32e0cSBarry Smith   ierr       = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr);
837d5b485abSSatish Balay   len        = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*);
83882502324SSatish Balay   ierr       = PetscMalloc(len,&sbuf2);CHKERRQ(ierr);
839d5b485abSSatish Balay   req_size   = (int*)(sbuf2 + nrqr);
840d5b485abSSatish Balay   req_source = req_size + nrqr;
841d5b485abSSatish Balay 
842d5b485abSSatish Balay   {
843df36731bSBarry Smith     Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
844ca161407SBarry Smith     int        *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
845d5b485abSSatish Balay 
846d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
847999d9058SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
848999d9058SBarry Smith       req_size[idex] = 0;
849999d9058SBarry Smith       rbuf1_i         = rbuf1[idex];
850d5b485abSSatish Balay       start           = 2*rbuf1_i[0] + 1;
851ca161407SBarry Smith       ierr            = MPI_Get_count(r_status1+i,MPI_INT,&end);CHKERRQ(ierr);
852999d9058SBarry Smith       ierr            = PetscMalloc(end*sizeof(int),&sbuf2[idex]);CHKERRQ(ierr);
853999d9058SBarry Smith       sbuf2_i         = sbuf2[idex];
854d5b485abSSatish Balay       for (j=start; j<end; j++) {
855d5b485abSSatish Balay         id               = rbuf1_i[j] - rstart;
856d5b485abSSatish Balay         ncols            = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
857d5b485abSSatish Balay         sbuf2_i[j]       = ncols;
858999d9058SBarry Smith         req_size[idex] += ncols;
859d5b485abSSatish Balay       }
860999d9058SBarry Smith       req_source[idex] = r_status1[i].MPI_SOURCE;
861d5b485abSSatish Balay       /* form the header */
862999d9058SBarry Smith       sbuf2_i[0]   = req_size[idex];
863d5b485abSSatish Balay       for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
864999d9058SBarry Smith       ierr = MPI_Isend(sbuf2_i,end,MPI_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr);
865d5b485abSSatish Balay     }
866d5b485abSSatish Balay   }
867606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
868606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
869d5b485abSSatish Balay 
870d5b485abSSatish Balay   /*  recv buffer sizes */
871d5b485abSSatish Balay   /* Receive messages*/
872d5b485abSSatish Balay 
873b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(int*),&rbuf3);CHKERRQ(ierr);
874b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr);
875b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr);
876b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr);
877b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr);
878d5b485abSSatish Balay 
879d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
880999d9058SBarry Smith     ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr);
881999d9058SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*sizeof(int),&rbuf3[idex]);CHKERRQ(ierr);
882999d9058SBarry Smith     ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr);
883999d9058SBarry Smith     ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPI_INT,
884999d9058SBarry Smith                      r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr);
885999d9058SBarry Smith     ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,
886999d9058SBarry Smith                      r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr);
887d5b485abSSatish Balay   }
888606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
889606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
890d5b485abSSatish Balay 
891d5b485abSSatish Balay   /* Wait on sends1 and sends2 */
892b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr);
893b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr);
894d5b485abSSatish Balay 
895ca161407SBarry Smith   ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);
896ca161407SBarry Smith   ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);
897606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
898606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
899606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
900606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
901d5b485abSSatish Balay 
902d5b485abSSatish Balay   /* Now allocate buffers for a->j, and send them off */
90382502324SSatish Balay   ierr = PetscMalloc((nrqr+1)*sizeof(int*),&sbuf_aj);CHKERRQ(ierr);
904d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
90582502324SSatish Balay   ierr = PetscMalloc((j+1)*sizeof(int),&sbuf_aj[0]);CHKERRQ(ierr);
906d5b485abSSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
907d5b485abSSatish Balay 
908b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr);
909d5b485abSSatish Balay   {
910d5b485abSSatish Balay      for (i=0; i<nrqr; i++) {
911d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
912d5b485abSSatish Balay       sbuf_aj_i = sbuf_aj[i];
913d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0] + 1;
914d5b485abSSatish Balay       ct2       = 0;
915d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
916d5b485abSSatish Balay         kmax = rbuf1[i][2*j];
917d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
918e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
919d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
920d5b485abSSatish Balay           ncols  = nzA + nzB;
921d5b485abSSatish Balay           cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
922d5b485abSSatish Balay 
923d5b485abSSatish Balay           /* load the column indices for this row into cols*/
924d5b485abSSatish Balay           cols  = sbuf_aj_i + ct2;
925d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
926d5b485abSSatish Balay             if ((ctmp = bmap[cworkB[l]]) < cstart)  cols[l] = ctmp;
927d5b485abSSatish Balay             else break;
928d5b485abSSatish Balay           }
929d5b485abSSatish Balay           imark = l;
930d5b485abSSatish Balay           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
931d5b485abSSatish Balay           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
932d5b485abSSatish Balay           ct2 += ncols;
933d5b485abSSatish Balay         }
934d5b485abSSatish Balay       }
935ca161407SBarry Smith       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
936d5b485abSSatish Balay     }
937d5b485abSSatish Balay   }
938b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr);
939b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr);
940d5b485abSSatish Balay 
941d5b485abSSatish Balay   /* Allocate buffers for a->a, and send them off */
94282502324SSatish Balay   ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr);
943d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
94482502324SSatish Balay   ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr);
945a2feddc5SSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
946d5b485abSSatish Balay 
947b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr);
948d5b485abSSatish Balay   {
949d5b485abSSatish Balay     for (i=0; i<nrqr; i++) {
950d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
951d5b485abSSatish Balay       sbuf_aa_i = sbuf_aa[i];
952d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0]+1;
953d5b485abSSatish Balay       ct2       = 0;
954d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
955d5b485abSSatish Balay         kmax = rbuf1_i[2*j];
956d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
957e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
958d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
959d5b485abSSatish Balay           ncols  = nzA + nzB;
960d5b485abSSatish Balay           cworkA = a_j + a_i[row];     cworkB = b_j + b_i[row];
961a2feddc5SSatish Balay           vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2;
962d5b485abSSatish Balay 
963d5b485abSSatish Balay           /* load the column values for this row into vals*/
9645b83ace0SSatish Balay           vals  = sbuf_aa_i+ct2*bs2;
965d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
9663a40ed3dSBarry Smith             if ((bmap[cworkB[l]]) < cstart) {
9673eda8832SBarry Smith               ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9683a40ed3dSBarry Smith             }
969d5b485abSSatish Balay             else break;
970d5b485abSSatish Balay           }
971d5b485abSSatish Balay           imark = l;
9723a40ed3dSBarry Smith           for (l=0; l<nzA; l++) {
9733eda8832SBarry Smith             ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9743a40ed3dSBarry Smith           }
9753a40ed3dSBarry Smith           for (l=imark; l<nzB; l++) {
9763eda8832SBarry Smith             ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
9773a40ed3dSBarry Smith           }
978d5b485abSSatish Balay           ct2 += ncols;
979d5b485abSSatish Balay         }
980d5b485abSSatish Balay       }
9813eda8832SBarry Smith       ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
982d5b485abSSatish Balay     }
983d5b485abSSatish Balay   }
984b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr);
985b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr);
986606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
987d5b485abSSatish Balay 
988d5b485abSSatish Balay   /* Form the matrix */
989d5b485abSSatish Balay   /* create col map */
990d5b485abSSatish Balay   {
991d5b485abSSatish Balay     int *icol_i;
992233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
993233c2ff6SSatish Balay     /* Create row map*/
994b0a32e0cSBarry Smith     ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&colmaps);CHKERRQ(ierr);
995ff0794f7SSatish Balay     for (i=0; i<ismax; i++) {
996ff0794f7SSatish Balay       ierr = PetscTableCreate(ncol[i]+1,&colmaps[i]);CHKERRQ(ierr);
997233c2ff6SSatish Balay     }
998233c2ff6SSatish Balay #else
999a2feddc5SSatish Balay     len     = (1+ismax)*sizeof(int*)+ ismax*c->Nbs*sizeof(int);
100082502324SSatish Balay     ierr    = PetscMalloc(len,&cmap);CHKERRQ(ierr);
1001d5b485abSSatish Balay     cmap[0] = (int*)(cmap + ismax);
1002549d3d68SSatish Balay     ierr    = PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(int));CHKERRQ(ierr);
1003a2feddc5SSatish Balay     for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; }
1004233c2ff6SSatish Balay #endif
1005d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1006d5b485abSSatish Balay       jmax   = ncol[i];
1007d5b485abSSatish Balay       icol_i = icol[i];
1008233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1009233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
1010233c2ff6SSatish Balay       for (j=0; j<jmax; j++) {
1011001aedefSBarry Smith 	ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr);
1012233c2ff6SSatish Balay       }
1013233c2ff6SSatish Balay #else
1014d5b485abSSatish Balay       cmap_i = cmap[i];
1015d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1016d5b485abSSatish Balay         cmap_i[icol_i[j]] = j+1;
1017d5b485abSSatish Balay       }
1018233c2ff6SSatish Balay #endif
1019d5b485abSSatish Balay     }
1020d5b485abSSatish Balay   }
1021d5b485abSSatish Balay 
1022d5b485abSSatish Balay   /* Create lens which is required for MatCreate... */
1023d5b485abSSatish Balay   for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
1024d5b485abSSatish Balay   len     = (1+ismax)*sizeof(int*)+ j*sizeof(int);
102582502324SSatish Balay   ierr    = PetscMalloc(len,&lens);CHKERRQ(ierr);
1026d5b485abSSatish Balay   lens[0] = (int*)(lens + ismax);
1027549d3d68SSatish Balay   ierr    = PetscMemzero(lens[0],j*sizeof(int));CHKERRQ(ierr);
1028d5b485abSSatish Balay   for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
1029d5b485abSSatish Balay 
1030d5b485abSSatish Balay   /* Update lens from local data */
1031d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1032d5b485abSSatish Balay     jmax   = nrow[i];
1033233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1034233c2ff6SSatish Balay     lcol1_gcol1 = colmaps[i];
1035233c2ff6SSatish Balay #else
1036d5b485abSSatish Balay     cmap_i = cmap[i];
1037233c2ff6SSatish Balay #endif
1038d5b485abSSatish Balay     irow_i = irow[i];
1039d5b485abSSatish Balay     lens_i = lens[i];
1040d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1041d5b485abSSatish Balay       row  = irow_i[j];
1042233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1043233c2ff6SSatish Balay       ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr);
1044233c2ff6SSatish Balay #else
1045d5b485abSSatish Balay       proc = rtable[row];
1046233c2ff6SSatish Balay #endif
1047d5b485abSSatish Balay       if (proc == rank) {
104836f4e84dSSatish Balay         /* Get indices from matA and then from matB */
104936f4e84dSSatish Balay         row    = row - rstart;
105036f4e84dSSatish Balay         nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
105136f4e84dSSatish Balay         cworkA =  a_j + a_i[row]; cworkB = b_j + b_i[row];
1052233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1053233c2ff6SSatish Balay        for (k=0; k<nzA; k++) {
1054233c2ff6SSatish Balay 	 ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1055233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1056233c2ff6SSatish Balay         }
1057233c2ff6SSatish Balay         for (k=0; k<nzB; k++) {
1058233c2ff6SSatish Balay 	  ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1059233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1060233c2ff6SSatish Balay         }
1061233c2ff6SSatish Balay #else
1062ca161407SBarry Smith         for (k=0; k<nzA; k++) {
106336f4e84dSSatish Balay           if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; }
1064ca161407SBarry Smith         }
1065ca161407SBarry Smith         for (k=0; k<nzB; k++) {
106636f4e84dSSatish Balay           if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; }
1067d5b485abSSatish Balay         }
1068233c2ff6SSatish Balay #endif
1069d5b485abSSatish Balay       }
1070a2feddc5SSatish Balay     }
1071ca161407SBarry Smith   }
1072233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1073d5b485abSSatish Balay   /* Create row map*/
1074b0a32e0cSBarry Smith   ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rowmaps);CHKERRQ(ierr);
1075ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
1076ff0794f7SSatish Balay     ierr = PetscTableCreate(nrow[i]+1,&rowmaps[i]);CHKERRQ(ierr);
1077233c2ff6SSatish Balay   }
1078233c2ff6SSatish Balay #else
1079233c2ff6SSatish Balay   /* Create row map*/
1080233c2ff6SSatish Balay   len     = (1+ismax)*sizeof(int*)+ ismax*Mbs*sizeof(int);
108182502324SSatish Balay   ierr    = PetscMalloc(len,&rmap);CHKERRQ(ierr);
1082d5b485abSSatish Balay   rmap[0] = (int*)(rmap + ismax);
1083233c2ff6SSatish Balay   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(int));CHKERRQ(ierr);
1084233c2ff6SSatish Balay   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;}
1085233c2ff6SSatish Balay #endif
1086d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1087d5b485abSSatish Balay     irow_i = irow[i];
1088d5b485abSSatish Balay     jmax   = nrow[i];
1089233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1090233c2ff6SSatish Balay     lrow1_grow1 = rowmaps[i];
1091233c2ff6SSatish Balay     for (j=0; j<jmax; j++) {
1092233c2ff6SSatish Balay       ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr);
1093233c2ff6SSatish Balay     }
1094233c2ff6SSatish Balay #else
1095233c2ff6SSatish Balay     rmap_i = rmap[i];
1096d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1097d5b485abSSatish Balay       rmap_i[irow_i[j]] = j;
1098d5b485abSSatish Balay     }
1099233c2ff6SSatish Balay #endif
1100d5b485abSSatish Balay   }
1101d5b485abSSatish Balay 
1102d5b485abSSatish Balay   /* Update lens from offproc data */
1103d5b485abSSatish Balay   {
1104d5b485abSSatish Balay     int *rbuf2_i,*rbuf3_i,*sbuf1_i;
1105d5b485abSSatish Balay 
1106d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1107ca161407SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits3,&i,r_status3+tmp2);CHKERRQ(ierr);
1108999d9058SBarry Smith       idex   = pa[i];
1109999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1110d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1111d5b485abSSatish Balay       ct1     = 2*jmax+1;
1112d5b485abSSatish Balay       ct2     = 0;
1113d5b485abSSatish Balay       rbuf2_i = rbuf2[i];
1114d5b485abSSatish Balay       rbuf3_i = rbuf3[i];
1115d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1116d5b485abSSatish Balay         is_no   = sbuf1_i[2*j-1];
1117d5b485abSSatish Balay         max1    = sbuf1_i[2*j];
1118d5b485abSSatish Balay         lens_i  = lens[is_no];
1119233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1120233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1121233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1122233c2ff6SSatish Balay #else
1123d5b485abSSatish Balay         cmap_i  = cmap[is_no];
1124d5b485abSSatish Balay 	rmap_i  = rmap[is_no];
1125233c2ff6SSatish Balay #endif
1126d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1127233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1128233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1129233c2ff6SSatish Balay           row--;
113029bbc08cSBarry Smith           if(row < 0) { SETERRQ(1,"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++) {
1136233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1137233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr);
1138233c2ff6SSatish Balay 	    if (tt) {
1139233c2ff6SSatish Balay               lens_i[row]++;
1140233c2ff6SSatish Balay             }
1141233c2ff6SSatish Balay #else
1142d5b485abSSatish Balay             if (cmap_i[rbuf3_i[ct2]]) {
1143d5b485abSSatish Balay               lens_i[row]++;
1144d5b485abSSatish Balay             }
1145233c2ff6SSatish Balay #endif
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);
1153ca161407SBarry Smith   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) {
1159d5b485abSSatish Balay     /*
1160d5b485abSSatish Balay         Assumes new rows are same length as the old rows, hence bug!
1161d5b485abSSatish Balay     */
1162d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1163df36731bSBarry Smith       mat = (Mat_SeqBAIJ *)(submats[i]->data);
116436f4e84dSSatish Balay       if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || mat->bs != bs)) {
116529bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1166d5b485abSSatish Balay       }
11670f5bd95cSBarry Smith       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(int),&flag);CHKERRQ(ierr);
11680f5bd95cSBarry Smith       if (flag == PETSC_FALSE) {
116929bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1170d5b485abSSatish Balay       }
1171d5b485abSSatish Balay       /* Initial matrix as if empty */
1172549d3d68SSatish Balay       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(int));CHKERRQ(ierr);
1173ce60de66SLois Curfman McInnes       submats[i]->factor = C->factor;
1174d5b485abSSatish Balay     }
1175ca161407SBarry Smith   } else {
1176d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1177e7ef3d9dSBarry Smith       ierr = MatCreate(PETSC_COMM_SELF,nrow[i]*bs,ncol[i]*bs,nrow[i]*bs,ncol[i]*bs,submats+i);CHKERRQ(ierr);
1178e7ef3d9dSBarry Smith       ierr = MatSetType(submats[i],A->type_name);CHKERRQ(ierr);
1179e7ef3d9dSBarry Smith       ierr = MatSeqBAIJSetPreallocation(submats[i],a->bs,0,lens[i]);CHKERRQ(ierr);
1180e7ef3d9dSBarry Smith       ierr = MatSeqSBAIJSetPreallocation(submats[i],a->bs,0,lens[i]);CHKERRQ(ierr);
1181d5b485abSSatish Balay     }
1182d5b485abSSatish Balay   }
1183d5b485abSSatish Balay 
1184d5b485abSSatish Balay   /* Assemble the matrices */
1185d5b485abSSatish Balay   /* First assemble the local rows */
1186d5b485abSSatish Balay   {
118736f4e84dSSatish Balay     int       ilen_row,*imat_ilen,*imat_j,*imat_i;
11883eda8832SBarry Smith     MatScalar *imat_a;
1189d5b485abSSatish Balay 
1190d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1191df36731bSBarry Smith       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1192d5b485abSSatish Balay       imat_ilen = mat->ilen;
1193d5b485abSSatish Balay       imat_j    = mat->j;
1194d5b485abSSatish Balay       imat_i    = mat->i;
1195d5b485abSSatish Balay       imat_a    = mat->a;
1196233c2ff6SSatish Balay 
1197233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1198233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
1199233c2ff6SSatish Balay       lrow1_grow1 = rowmaps[i];
1200233c2ff6SSatish Balay #else
1201d5b485abSSatish Balay       cmap_i    = cmap[i];
1202d5b485abSSatish Balay       rmap_i    = rmap[i];
1203233c2ff6SSatish Balay #endif
1204d5b485abSSatish Balay       irow_i    = irow[i];
1205d5b485abSSatish Balay       jmax      = nrow[i];
1206d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1207d5b485abSSatish Balay         row      = irow_i[j];
1208233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1209233c2ff6SSatish Balay 	ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr);
1210233c2ff6SSatish Balay #else
1211d5b485abSSatish Balay 	proc = rtable[row];
1212233c2ff6SSatish Balay #endif
1213d5b485abSSatish Balay         if (proc == rank) {
121436f4e84dSSatish Balay           row      = row - rstart;
121536f4e84dSSatish Balay           nzA      = a_i[row+1] - a_i[row];
121636f4e84dSSatish Balay           nzB      = b_i[row+1] - b_i[row];
121736f4e84dSSatish Balay           cworkA   = a_j + a_i[row];
121836f4e84dSSatish Balay           cworkB   = b_j + b_i[row];
121936f4e84dSSatish Balay           vworkA   = a_a + a_i[row]*bs2;
122036f4e84dSSatish Balay           vworkB   = b_a + b_i[row]*bs2;
1221233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1222233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr);
1223233c2ff6SSatish Balay           row--;
122429bbc08cSBarry Smith           if (row < 0) { SETERRQ(1,"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];
122936f4e84dSSatish Balay           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*/
123436f4e84dSSatish Balay           for (l=0; l<nzB; l++) {
123536f4e84dSSatish Balay 	    if ((ctmp = bmap[cworkB[l]]) < cstart) {
1236233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1237233c2ff6SSatish Balay 	      ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr);
1238233c2ff6SSatish Balay 	      if (tcol) {
1239233c2ff6SSatish Balay #else
124036f4e84dSSatish Balay               if ((tcol = cmap_i[ctmp])) {
1241233c2ff6SSatish Balay #endif
1242df36731bSBarry Smith                 *mat_j++ = tcol - 1;
12433eda8832SBarry Smith                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1244549d3d68SSatish Balay                 mat_a   += bs2;
1245d5b485abSSatish Balay                 ilen_row++;
1246d5b485abSSatish Balay               }
1247ca161407SBarry Smith             } else break;
124836f4e84dSSatish Balay           }
124936f4e84dSSatish Balay           imark = l;
125036f4e84dSSatish Balay           for (l=0; l<nzA; l++) {
1251233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1252233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1253233c2ff6SSatish Balay 	    if (tcol) {
1254233c2ff6SSatish Balay #else
125536f4e84dSSatish Balay 	    if ((tcol = cmap_i[cstart + cworkA[l]])) {
1256233c2ff6SSatish Balay #endif
125736f4e84dSSatish Balay               *mat_j++ = tcol - 1;
12583eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1259549d3d68SSatish Balay               mat_a   += bs2;
126036f4e84dSSatish Balay               ilen_row++;
126136f4e84dSSatish Balay             }
126236f4e84dSSatish Balay           }
126336f4e84dSSatish Balay           for (l=imark; l<nzB; l++) {
1264233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1265233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1266233c2ff6SSatish Balay 	    if (tcol) {
1267233c2ff6SSatish Balay #else
126836f4e84dSSatish Balay             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1269233c2ff6SSatish Balay #endif
127036f4e84dSSatish Balay               *mat_j++ = tcol - 1;
12713eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1272549d3d68SSatish Balay               mat_a   += bs2;
127336f4e84dSSatish Balay               ilen_row++;
127436f4e84dSSatish Balay             }
127536f4e84dSSatish Balay           }
1276d5b485abSSatish Balay           imat_ilen[row] = ilen_row;
1277d5b485abSSatish Balay         }
1278d5b485abSSatish Balay       }
127936f4e84dSSatish Balay 
1280d5b485abSSatish Balay     }
1281d5b485abSSatish Balay   }
1282d5b485abSSatish Balay 
1283d5b485abSSatish Balay   /*   Now assemble the off proc rows*/
1284d5b485abSSatish Balay   {
1285d5b485abSSatish Balay     int       *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1286d5b485abSSatish Balay     int       *imat_j,*imat_i;
12873eda8832SBarry Smith     MatScalar *imat_a,*rbuf4_i;
1288d5b485abSSatish Balay 
1289d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1290ca161407SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits4,&i,r_status4+tmp2);CHKERRQ(ierr);
1291999d9058SBarry Smith       idex   = pa[i];
1292999d9058SBarry Smith       sbuf1_i = sbuf1[idex];
1293d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1294d5b485abSSatish Balay       ct1     = 2*jmax + 1;
1295d5b485abSSatish Balay       ct2     = 0;
1296d5b485abSSatish Balay       rbuf2_i = rbuf2[i];
1297d5b485abSSatish Balay       rbuf3_i = rbuf3[i];
1298d5b485abSSatish Balay       rbuf4_i = rbuf4[i];
1299d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1300d5b485abSSatish Balay         is_no     = sbuf1_i[2*j-1];
1301233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1302233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1303233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1304233c2ff6SSatish Balay #else
1305d5b485abSSatish Balay         rmap_i    = rmap[is_no];
1306d5b485abSSatish Balay         cmap_i    = cmap[is_no];
1307233c2ff6SSatish Balay #endif
1308df36731bSBarry Smith         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1309d5b485abSSatish Balay         imat_ilen = mat->ilen;
1310d5b485abSSatish Balay         imat_j    = mat->j;
1311d5b485abSSatish Balay         imat_i    = mat->i;
1312d5b485abSSatish Balay         imat_a    = mat->a;
1313d5b485abSSatish Balay         max1      = sbuf1_i[2*j];
1314d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1315d5b485abSSatish Balay           row   = sbuf1_i[ct1];
1316233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1317233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr);
1318233c2ff6SSatish Balay           row--;
131929bbc08cSBarry Smith           if(row < 0) { SETERRQ(1,"row not found in table"); }
1320233c2ff6SSatish Balay #else
1321d5b485abSSatish Balay           row   = rmap_i[row];
1322233c2ff6SSatish Balay #endif
1323d5b485abSSatish Balay           ilen  = imat_ilen[row];
1324df36731bSBarry Smith           mat_i = imat_i[row];
132536f4e84dSSatish Balay           mat_a = imat_a + mat_i*bs2;
1326d5b485abSSatish Balay           mat_j = imat_j + mat_i;
1327d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1328d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1329233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1330233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1331233c2ff6SSatish Balay 	    if (tcol) {
1332233c2ff6SSatish Balay #else
1333d5b485abSSatish Balay 	    if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1334233c2ff6SSatish Balay #endif
1335df36731bSBarry Smith               *mat_j++    = tcol - 1;
133636f4e84dSSatish Balay               /* *mat_a++= rbuf4_i[ct2]; */
13373eda8832SBarry Smith               ierr        = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1338549d3d68SSatish Balay               mat_a      += bs2;
1339d5b485abSSatish Balay               ilen++;
1340d5b485abSSatish Balay             }
1341d5b485abSSatish Balay           }
1342d5b485abSSatish Balay           imat_ilen[row] = ilen;
1343d5b485abSSatish Balay         }
1344d5b485abSSatish Balay       }
1345d5b485abSSatish Balay     }
1346d5b485abSSatish Balay   }
1347606d414cSSatish Balay   ierr = PetscFree(r_status4);CHKERRQ(ierr);
1348606d414cSSatish Balay   ierr = PetscFree(r_waits4);CHKERRQ(ierr);
1349ca161407SBarry Smith   ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);
1350606d414cSSatish Balay   ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1351606d414cSSatish Balay   ierr = PetscFree(s_status4);CHKERRQ(ierr);
1352d5b485abSSatish Balay 
1353d5b485abSSatish Balay   /* Restore the indices */
1354d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1355d5b485abSSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1356d5b485abSSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1357d5b485abSSatish Balay   }
1358d5b485abSSatish Balay 
1359d5b485abSSatish Balay   /* Destroy allocated memory */
1360606d414cSSatish Balay   ierr = PetscFree(irow);CHKERRQ(ierr);
1361606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
1362606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
1363d5b485abSSatish Balay 
1364606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1365606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1366d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
1367606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1368d5b485abSSatish Balay   }
1369d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1370606d414cSSatish Balay     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1371606d414cSSatish Balay     ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);
1372d5b485abSSatish Balay   }
1373d5b485abSSatish Balay 
1374606d414cSSatish Balay   ierr = PetscFree(sbuf2);CHKERRQ(ierr);
1375606d414cSSatish Balay   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1376606d414cSSatish Balay   ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1377606d414cSSatish Balay   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1378606d414cSSatish Balay   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
1379606d414cSSatish Balay   ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1380606d414cSSatish Balay   ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
1381d5b485abSSatish Balay 
1382233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1383ff0794f7SSatish Balay   for (i=0; i<ismax; i++){
1384233c2ff6SSatish Balay     ierr = PetscTableDelete(rowmaps[i]);CHKERRQ(ierr);
1385233c2ff6SSatish Balay     ierr = PetscTableDelete(colmaps[i]);CHKERRQ(ierr);
1386233c2ff6SSatish Balay   }
1387233c2ff6SSatish Balay   ierr = PetscFree(colmaps);CHKERRQ(ierr);
1388233c2ff6SSatish Balay   ierr = PetscFree(rowmaps);CHKERRQ(ierr);
1389233c2ff6SSatish Balay #else
1390606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
1391233c2ff6SSatish Balay   ierr = PetscFree(cmap);CHKERRQ(ierr);
1392233c2ff6SSatish Balay #endif
1393606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1394d5b485abSSatish Balay 
1395d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
139636f4e84dSSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139736f4e84dSSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1398d5b485abSSatish Balay   }
139934e6ae18SSatish Balay 
14003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1401d5b485abSSatish Balay }
1402ca161407SBarry Smith 
1403