xref: /petsc/src/mat/impls/baij/mpi/baijov.c (revision b0a32e0c6855ee6a6cd3495fa7da12ea9885bc5d)
1*b0a32e0cSBarry Smith /*$Id: baijov.c,v 1.55 2000/09/28 21:11:32 bsmith Exp bsmith $*/
2ca161407SBarry Smith 
3d5b485abSSatish Balay /*
4d5b485abSSatish Balay    Routines to compute overlapping regions of a parallel MPI matrix
5d5b485abSSatish Balay   and to find submatrices that were shared across processors.
6d5b485abSSatish Balay */
770f55243SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h"
80a835dfdSSatish Balay #include "petscbt.h"
9d5b485abSSatish Balay 
10df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat,int,IS *);
11df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Local(Mat,int,char **,int*,int**);
12df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat,int,int **,int**,int*);
13ca44d042SBarry Smith EXTERN int MatGetRow_MPIBAIJ(Mat,int,int*,int**,Scalar**);
14ca44d042SBarry Smith EXTERN int MatRestoreRow_MPIBAIJ(Mat,int,int*,int**,Scalar**);
15d5b485abSSatish Balay 
165615d1e5SSatish Balay #undef __FUNC__
17*b0a32e0cSBarry Smith #define __FUNC__ "MatCompressIndicesGeneral_MPIBAIJ"
18e757655aSSatish Balay static int MatCompressIndicesGeneral_MPIBAIJ(Mat C,int imax,IS *is_in,IS *is_out)
19df36731bSBarry Smith {
20df36731bSBarry Smith   Mat_MPIBAIJ  *baij = (Mat_MPIBAIJ*)C->data;
21233c2ff6SSatish Balay   int          ierr,isz,bs = baij->bs,n,i,j,*idx,ival;
22233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
23233c2ff6SSatish Balay   PetscTable   gid1_lid1;
24233c2ff6SSatish Balay   int tt, gid1, *nidx;
25233c2ff6SSatish Balay   PetscTablePosition tpos;
26233c2ff6SSatish Balay #else
27233c2ff6SSatish Balay   int          Nbs,*nidx;
28f1af5d2fSBarry Smith   PetscBT      table;
29233c2ff6SSatish Balay #endif
30df36731bSBarry Smith 
313a40ed3dSBarry Smith   PetscFunctionBegin;
32233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
33233c2ff6SSatish Balay   ierr = PetscTableCreate(baij->mbs,&gid1_lid1);CHKERRQ(ierr);
34233c2ff6SSatish Balay #else
35df36731bSBarry Smith   Nbs   = baij->Nbs;
36*b0a32e0cSBarry Smith ierr = PetscMalloc((Nbs+1)*sizeof(int),&  nidx  );CHKERRQ(ierr);
376831982aSBarry Smith   ierr  = PetscBTCreate(Nbs,table);CHKERRQ(ierr);
38233c2ff6SSatish Balay #endif
39df36731bSBarry Smith   for (i=0; i<imax; i++) {
40df36731bSBarry Smith     isz  = 0;
41233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
42233c2ff6SSatish Balay     ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
43233c2ff6SSatish Balay #else
446831982aSBarry Smith     ierr = PetscBTMemzero(Nbs,table);CHKERRQ(ierr);
45233c2ff6SSatish Balay #endif
4636f4e84dSSatish Balay     ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr);
47b9b97703SBarry Smith     ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr);
48e757655aSSatish Balay     for (j=0; j<n ; j++) {
49df36731bSBarry Smith       ival = idx[j]/bs; /* convert the indices into block indices */
50233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
51233c2ff6SSatish Balay       ierr = PetscTableFind(gid1_lid1,ival+1,&tt);CHKERRQ(ierr);
52233c2ff6SSatish Balay       if (!tt) {
53233c2ff6SSatish Balay 	ierr = PetscTableAdd(gid1_lid1,ival+1,isz+1);CHKERRQ(ierr);
54233c2ff6SSatish Balay         isz++;
55233c2ff6SSatish Balay       }
56233c2ff6SSatish Balay #else
5729bbc08cSBarry Smith       if (ival>Nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
58f1af5d2fSBarry Smith       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
59233c2ff6SSatish Balay #endif
60df36731bSBarry Smith     }
6136f4e84dSSatish Balay     ierr = ISRestoreIndices(is_in[i],&idx);CHKERRQ(ierr);
62233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
63*b0a32e0cSBarry Smith ierr = PetscMalloc((isz+1)*sizeof(int),&    nidx );CHKERRQ(ierr);
64233c2ff6SSatish Balay     ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
65233c2ff6SSatish Balay     j = 0;
66233c2ff6SSatish Balay     while (tpos) {
67233c2ff6SSatish Balay       ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid1,&tt);CHKERRQ(ierr);
6829bbc08cSBarry Smith       if (tt-- > isz) { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than array-dim"); }
69233c2ff6SSatish Balay       nidx[tt] = gid1 - 1;
70233c2ff6SSatish Balay       j++;
71df36731bSBarry Smith     }
7229bbc08cSBarry Smith     if (j != isz) { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"table error: jj != isz"); }
73233c2ff6SSatish Balay     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is_out+i));CHKERRQ(ierr);
74233c2ff6SSatish Balay     ierr = PetscFree(nidx);CHKERRQ(ierr);
75233c2ff6SSatish Balay #else
76233c2ff6SSatish Balay     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is_out+i));CHKERRQ(ierr);
77233c2ff6SSatish Balay #endif
78233c2ff6SSatish Balay   }
79233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
80233c2ff6SSatish Balay   ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr);
81233c2ff6SSatish Balay #else
826831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
83606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
84233c2ff6SSatish Balay #endif
853a40ed3dSBarry Smith   PetscFunctionReturn(0);
86df36731bSBarry Smith }
87df36731bSBarry Smith 
885615d1e5SSatish Balay #undef __FUNC__
89*b0a32e0cSBarry Smith #define __FUNC__ "MatCompressIndicesSorted_MPIBAIJ"
90e757655aSSatish Balay static int MatCompressIndicesSorted_MPIBAIJ(Mat C,int imax,IS *is_in,IS *is_out)
91e757655aSSatish Balay {
92e757655aSSatish Balay   Mat_MPIBAIJ  *baij = (Mat_MPIBAIJ*)C->data;
93233c2ff6SSatish Balay   int          ierr,bs=baij->bs,i,j,k,val,n,*idx,*nidx,*idx_local;
94e757655aSSatish Balay   PetscTruth   flg;
95233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
96233c2ff6SSatish Balay   int maxsz;
97233c2ff6SSatish Balay #else
98233c2ff6SSatish Balay   int Nbs=baij->Nbs;
99233c2ff6SSatish Balay #endif
1003a40ed3dSBarry Smith   PetscFunctionBegin;
101e757655aSSatish Balay   for (i=0; i<imax; i++) {
102e757655aSSatish Balay     ierr = ISSorted(is_in[i],&flg);CHKERRQ(ierr);
10329bbc08cSBarry Smith     if (!flg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Indices are not sorted");
104e757655aSSatish Balay   }
105233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
106233c2ff6SSatish Balay   /* Now check max size */
107233c2ff6SSatish Balay   for (i=0,maxsz=0; i<imax; i++) {
108233c2ff6SSatish Balay     ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr);
109233c2ff6SSatish Balay     ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr);
11029bbc08cSBarry Smith     if (n%bs !=0) SETERRA(1,"Indices are not block ordered");
111233c2ff6SSatish Balay     n = n/bs; /* The reduced index size */
112233c2ff6SSatish Balay     if (n > maxsz) maxsz = n;
113233c2ff6SSatish Balay   }
114*b0a32e0cSBarry Smith ierr = PetscMalloc((maxsz+1)*sizeof(int),&  nidx  );CHKERRQ(ierr);
115233c2ff6SSatish Balay #else
116*b0a32e0cSBarry Smith ierr = PetscMalloc((Nbs+1)*sizeof(int),&  nidx  );CHKERRQ(ierr);
117233c2ff6SSatish Balay #endif
1183a40ed3dSBarry Smith   /* Now check if the indices are in block order */
119e757655aSSatish Balay   for (i=0; i<imax; i++) {
120e757655aSSatish Balay     ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr);
121b9b97703SBarry Smith     ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr);
12229bbc08cSBarry Smith     if (n%bs !=0) SETERRA(1,"Indices are not block ordered");
123e757655aSSatish Balay 
124e757655aSSatish Balay     n = n/bs; /* The reduced index size */
125e757655aSSatish Balay     idx_local = idx;
126e757655aSSatish Balay     for (j=0; j<n ; j++) {
127e757655aSSatish Balay       val = idx_local[0];
12829bbc08cSBarry Smith       if (val%bs != 0) SETERRA(1,"Indices are not block ordered");
129e757655aSSatish Balay       for (k=0; k<bs; k++) {
13029bbc08cSBarry Smith         if (val+k != idx_local[k]) SETERRA(1,"Indices are not block ordered");
131e757655aSSatish Balay       }
132e757655aSSatish Balay       nidx[j] = val/bs;
133e757655aSSatish Balay       idx_local +=bs;
134e757655aSSatish Balay     }
135e757655aSSatish Balay     ierr = ISRestoreIndices(is_in[i],&idx);CHKERRQ(ierr);
136029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,n,nidx,(is_out+i));CHKERRQ(ierr);
137e757655aSSatish Balay   }
138606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
139233c2ff6SSatish Balay 
1403a40ed3dSBarry Smith   PetscFunctionReturn(0);
141e757655aSSatish Balay }
142e757655aSSatish Balay 
1435615d1e5SSatish Balay #undef __FUNC__
144*b0a32e0cSBarry Smith #define __FUNC__ "MatExpandIndices_MPIBAIJ"
14536f4e84dSSatish Balay static int MatExpandIndices_MPIBAIJ(Mat C,int imax,IS *is_in,IS *is_out)
146df36731bSBarry Smith {
147df36731bSBarry Smith   Mat_MPIBAIJ  *baij = (Mat_MPIBAIJ*)C->data;
148233c2ff6SSatish Balay   int          ierr,bs = baij->bs,n,i,j,k,*idx,*nidx;
149233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
150233c2ff6SSatish Balay   int maxsz;
151233c2ff6SSatish Balay #else
152233c2ff6SSatish Balay   int Nbs = baij->Nbs;
153233c2ff6SSatish Balay #endif
1543a40ed3dSBarry Smith   PetscFunctionBegin;
155df36731bSBarry Smith 
156233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
157233c2ff6SSatish Balay   /* Now check max size */
158233c2ff6SSatish Balay   for (i=0,maxsz=0; i<imax; i++) {
159233c2ff6SSatish Balay     ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr);
160233c2ff6SSatish Balay     ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr);
161233c2ff6SSatish Balay     if (n*bs > maxsz) maxsz = n*bs;
162233c2ff6SSatish Balay   }
163*b0a32e0cSBarry Smith ierr = PetscMalloc((maxsz+1)*sizeof(int),&  nidx  );CHKERRQ(ierr);
164233c2ff6SSatish Balay #else
165*b0a32e0cSBarry Smith ierr = PetscMalloc((Nbs*bs+1)*sizeof(int),&  nidx  );CHKERRQ(ierr);
166233c2ff6SSatish Balay #endif
167df36731bSBarry Smith 
168df36731bSBarry Smith   for (i=0; i<imax; i++) {
16936f4e84dSSatish Balay     ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr);
170b9b97703SBarry Smith     ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr);
171df36731bSBarry Smith     for (j=0; j<n ; ++j){
172df36731bSBarry Smith       for (k=0; k<bs; k++)
173df36731bSBarry Smith         nidx[j*bs+k] = idx[j]*bs+k;
174df36731bSBarry Smith     }
17536f4e84dSSatish Balay     ierr = ISRestoreIndices(is_in[i],&idx);CHKERRQ(ierr);
176329f5518SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,n*bs,nidx,is_out+i);CHKERRQ(ierr);
177df36731bSBarry Smith   }
178606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
1793a40ed3dSBarry Smith   PetscFunctionReturn(0);
180df36731bSBarry Smith }
181df36731bSBarry Smith 
182df36731bSBarry Smith 
1835615d1e5SSatish Balay #undef __FUNC__
184*b0a32e0cSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ"
185df36731bSBarry Smith int MatIncreaseOverlap_MPIBAIJ(Mat C,int imax,IS *is,int ov)
186d5b485abSSatish Balay {
187d5b485abSSatish Balay   int i,ierr;
18836f4e84dSSatish Balay   IS  *is_new;
18936f4e84dSSatish Balay 
1903a40ed3dSBarry Smith   PetscFunctionBegin;
191*b0a32e0cSBarry Smith ierr = PetscMalloc(imax*sizeof(IS),&(  is_new ));CHKERRQ(ierr);
192df36731bSBarry Smith   /* Convert the indices into block format */
193e757655aSSatish Balay   ierr = MatCompressIndicesGeneral_MPIBAIJ(C,imax,is,is_new);CHKERRQ(ierr);
19429bbc08cSBarry Smith   if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
195d5b485abSSatish Balay   for (i=0; i<ov; ++i) {
19636f4e84dSSatish Balay     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr);
197d5b485abSSatish Balay   }
198ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);}
19936f4e84dSSatish Balay   ierr = MatExpandIndices_MPIBAIJ(C,imax,is_new,is);CHKERRQ(ierr);
200ca161407SBarry Smith   for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);}
201606d414cSSatish Balay   ierr = PetscFree(is_new);CHKERRQ(ierr);
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
203d5b485abSSatish Balay }
204d5b485abSSatish Balay 
205d5b485abSSatish Balay /*
206d5b485abSSatish Balay   Sample message format:
207d5b485abSSatish Balay   If a processor A wants processor B to process some elements corresponding
208d5b485abSSatish Balay   to index sets 1s[1], is[5]
209d5b485abSSatish Balay   mesg [0] = 2   (no of index sets in the mesg)
210d5b485abSSatish Balay   -----------
211d5b485abSSatish Balay   mesg [1] = 1 => is[1]
212d5b485abSSatish Balay   mesg [2] = sizeof(is[1]);
213d5b485abSSatish Balay   -----------
214d5b485abSSatish Balay   mesg [5] = 5  => is[5]
215d5b485abSSatish Balay   mesg [6] = sizeof(is[5]);
216d5b485abSSatish Balay   -----------
217d5b485abSSatish Balay   mesg [7]
218d5b485abSSatish Balay   mesg [n]  datas[1]
219d5b485abSSatish Balay   -----------
220d5b485abSSatish Balay   mesg[n+1]
221d5b485abSSatish Balay   mesg[m]  data(is[5])
222d5b485abSSatish Balay   -----------
223d5b485abSSatish Balay 
224d5b485abSSatish Balay   Notes:
225d5b485abSSatish Balay   nrqs - no of requests sent (or to be sent out)
226d5b485abSSatish Balay   nrqr - no of requests recieved (which have to be or which have been processed
227d5b485abSSatish Balay */
2285615d1e5SSatish Balay #undef __FUNC__
229*b0a32e0cSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Once"
230df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat C,int imax,IS *is)
231d5b485abSSatish Balay {
232df36731bSBarry Smith   Mat_MPIBAIJ  *c = (Mat_MPIBAIJ*)C->data;
233d5b485abSSatish Balay   int         **idx,*n,*w1,*w2,*w3,*w4,*rtable,**data,len,*idx_i;
234df36731bSBarry Smith   int         size,rank,Mbs,i,j,k,ierr,**rbuf,row,proc,nrqs,msz,**outdat,**ptr;
235d5b485abSSatish Balay   int         *ctr,*pa,tag,*tmp,bsz,nrqr,*isz,*isz1,**xdata,bsz1,**rbuf2;
236f1af5d2fSBarry Smith   PetscBT     *table;
237d5b485abSSatish Balay   MPI_Comm    comm;
238d5b485abSSatish Balay   MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2;
239d5b485abSSatish Balay   MPI_Status  *s_status,*recv_status;
240d5b485abSSatish Balay 
2413a40ed3dSBarry Smith   PetscFunctionBegin;
242d5b485abSSatish Balay   comm   = C->comm;
243d5b485abSSatish Balay   tag    = C->tag;
244d5b485abSSatish Balay   size   = c->size;
245d5b485abSSatish Balay   rank   = c->rank;
246df36731bSBarry Smith   Mbs    = c->Mbs;
247d5b485abSSatish Balay 
248df36731bSBarry Smith   len    = (imax+1)*sizeof(int*)+ (imax + Mbs)*sizeof(int);
249*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  idx    ));CHKERRQ(ierr);
250d5b485abSSatish Balay   n      = (int*)(idx + imax);
251d5b485abSSatish Balay   rtable = n + imax;
252d5b485abSSatish Balay 
253d5b485abSSatish Balay   for (i=0; i<imax; i++) {
254d5b485abSSatish Balay     ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr);
255b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
256d5b485abSSatish Balay   }
257d5b485abSSatish Balay 
258d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
259d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
260d5b485abSSatish Balay     len = c->rowners[i+1];
261d5b485abSSatish Balay     for (; j<len; j++) {
262d5b485abSSatish Balay       rtable[j] = i;
263d5b485abSSatish Balay     }
264d5b485abSSatish Balay   }
265d5b485abSSatish Balay 
266d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
267d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
268*b0a32e0cSBarry Smith ierr = PetscMalloc(size*4*sizeof(int),&(  w1   ));CHKERRQ(ierr);/*  mesg size */
269d5b485abSSatish Balay   w2   = w1 + size;       /* if w2[i] marked, then a message to proc i*/
270d5b485abSSatish Balay   w3   = w2 + size;       /* no of IS that needs to be sent to proc i */
271d5b485abSSatish Balay   w4   = w3 + size;       /* temp work space used in determining w1, w2, w3 */
272549d3d68SSatish Balay   ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/
273d5b485abSSatish Balay   for (i=0; i<imax; i++) {
274549d3d68SSatish Balay     ierr  = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/
275d5b485abSSatish Balay     idx_i = idx[i];
276d5b485abSSatish Balay     len   = n[i];
277d5b485abSSatish Balay     for (j=0; j<len; j++) {
278d5b485abSSatish Balay       row  = idx_i[j];
2796b41c64dSBarry Smith       if (row < 0) {
28029bbc08cSBarry Smith         SETERRQ(1,"Index set cannot have negative entries");
2816b41c64dSBarry Smith       }
282d5b485abSSatish Balay       proc = rtable[row];
283d5b485abSSatish Balay       w4[proc]++;
284d5b485abSSatish Balay     }
285d5b485abSSatish Balay     for (j=0; j<size; j++){
286d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
287d5b485abSSatish Balay     }
288d5b485abSSatish Balay   }
289d5b485abSSatish Balay 
290d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
291d5b485abSSatish Balay   msz      = 0;              /* total mesg length (for all proc */
292d5b485abSSatish Balay   w1[rank] = 0;              /* no mesg sent to intself */
293d5b485abSSatish Balay   w3[rank] = 0;
294d5b485abSSatish Balay   for (i=0; i<size; i++) {
295d5b485abSSatish Balay     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
296d5b485abSSatish Balay   }
297d5b485abSSatish Balay   /* pa - is list of processors to communicate with */
298*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int),&  pa );CHKERRQ(ierr);
299d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
300d5b485abSSatish Balay     if (w1[i]) {pa[j] = i; j++;}
301d5b485abSSatish Balay   }
302d5b485abSSatish Balay 
303d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
304d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
305d5b485abSSatish Balay     j      = pa[i];
306d5b485abSSatish Balay     w1[j] += w2[j] + 2*w3[j];
307d5b485abSSatish Balay     msz   += w1[j];
308d5b485abSSatish Balay   }
309d5b485abSSatish Balay 
310d5b485abSSatish Balay 
311d5b485abSSatish Balay   /* Do a global reduction to determine how many messages to expect*/
312d5b485abSSatish Balay   {
3136831982aSBarry Smith     int *rw1;
314*b0a32e0cSBarry Smith ierr = PetscMalloc(2*size*sizeof(int),&(    rw1   ));CHKERRQ(ierr);
3156831982aSBarry Smith     ierr  = MPI_Allreduce(w1,rw1,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr);
316d5b485abSSatish Balay     bsz   = rw1[rank];
3176831982aSBarry Smith     nrqr  = rw1[size+rank];
318606d414cSSatish Balay     ierr  = PetscFree(rw1);CHKERRQ(ierr);
319d5b485abSSatish Balay   }
320d5b485abSSatish Balay 
321d5b485abSSatish Balay   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
322d5b485abSSatish Balay   len     = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
323*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  rbuf    ));CHKERRQ(ierr);
324d5b485abSSatish Balay   rbuf[0] = (int*)(rbuf + nrqr);
325d5b485abSSatish Balay   for (i=1; i<nrqr; ++i) rbuf[i] = rbuf[i-1] + bsz;
326d5b485abSSatish Balay 
327d5b485abSSatish Balay   /* Post the receives */
328*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&  r_waits1 );CHKERRQ(ierr);
329d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
330ca161407SBarry Smith     ierr = MPI_Irecv(rbuf[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits1+i);CHKERRQ(ierr);
331d5b485abSSatish Balay   }
332d5b485abSSatish Balay 
333d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
334d5b485abSSatish Balay   len    = 2*size*sizeof(int*) + (size+msz)*sizeof(int);
335*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  outdat ));CHKERRQ(ierr);
336d5b485abSSatish Balay   ptr    = outdat + size;     /* Pointers to the data in outgoing buffers */
337549d3d68SSatish Balay   ierr   = PetscMemzero(outdat,2*size*sizeof(int*));CHKERRQ(ierr);
338d5b485abSSatish Balay   tmp    = (int*)(outdat + 2*size);
339d5b485abSSatish Balay   ctr    = tmp + msz;
340d5b485abSSatish Balay 
341d5b485abSSatish Balay   {
342d5b485abSSatish Balay     int *iptr = tmp,ict  = 0;
343d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
344d5b485abSSatish Balay       j         = pa[i];
345d5b485abSSatish Balay       iptr     +=  ict;
346d5b485abSSatish Balay       outdat[j] = iptr;
347d5b485abSSatish Balay       ict       = w1[j];
348d5b485abSSatish Balay     }
349d5b485abSSatish Balay   }
350d5b485abSSatish Balay 
351d5b485abSSatish Balay   /* Form the outgoing messages */
352d5b485abSSatish Balay   /*plug in the headers*/
353d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
354d5b485abSSatish Balay     j            = pa[i];
355d5b485abSSatish Balay     outdat[j][0] = 0;
356549d3d68SSatish Balay     ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr);
357d5b485abSSatish Balay     ptr[j]       = outdat[j] + 2*w3[j] + 1;
358d5b485abSSatish Balay   }
359d5b485abSSatish Balay 
360d5b485abSSatish Balay   /* Memory for doing local proc's work*/
361d5b485abSSatish Balay   {
362d5b485abSSatish Balay     int  *d_p;
363d5b485abSSatish Balay     char *t_p;
364d5b485abSSatish Balay 
365f1af5d2fSBarry Smith     len   = (imax)*(sizeof(PetscBT) + sizeof(int*)+ sizeof(int)) +
3662a8f2162SSatish Balay       (Mbs)*imax*sizeof(int)  + (Mbs/BITSPERBYTE+1)*imax*sizeof(char) + 1;
367*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(    table ));CHKERRQ(ierr);
368549d3d68SSatish Balay     ierr  = PetscMemzero(table,len);CHKERRQ(ierr);
369d5b485abSSatish Balay     data  = (int **)(table + imax);
370bbd702dbSSatish Balay     isz   = (int  *)(data  + imax);
371bbd702dbSSatish Balay     d_p   = (int  *)(isz   + imax);
372bbd702dbSSatish Balay     t_p   = (char *)(d_p   + Mbs*imax);
373bbd702dbSSatish Balay     for (i=0; i<imax; i++) {
3742a8f2162SSatish Balay       table[i] = t_p + (Mbs/BITSPERBYTE+1)*i;
375bbd702dbSSatish Balay       data[i]  = d_p + (Mbs)*i;
376d5b485abSSatish Balay     }
377d5b485abSSatish Balay   }
378d5b485abSSatish Balay 
379d5b485abSSatish Balay   /* Parse the IS and update local tables and the outgoing buf with the data*/
380d5b485abSSatish Balay   {
381d5b485abSSatish Balay     int     n_i,*data_i,isz_i,*outdat_j,ctr_j;
382f1af5d2fSBarry Smith     PetscBT table_i;
383d5b485abSSatish Balay 
384d5b485abSSatish Balay     for (i=0; i<imax; i++) {
385549d3d68SSatish Balay       ierr    = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr);
386d5b485abSSatish Balay       n_i     = n[i];
387d5b485abSSatish Balay       table_i = table[i];
388d5b485abSSatish Balay       idx_i   = idx[i];
389d5b485abSSatish Balay       data_i  = data[i];
390d5b485abSSatish Balay       isz_i   = isz[i];
391d5b485abSSatish Balay       for (j=0;  j<n_i; j++) {  /* parse the indices of each IS */
392d5b485abSSatish Balay         row  = idx_i[j];
393d5b485abSSatish Balay         proc = rtable[row];
394d5b485abSSatish Balay         if (proc != rank) { /* copy to the outgoing buffer */
395d5b485abSSatish Balay           ctr[proc]++;
396d5b485abSSatish Balay           *ptr[proc] = row;
397d5b485abSSatish Balay           ptr[proc]++;
398d5b485abSSatish Balay         }
399d5b485abSSatish Balay         else { /* Update the local table */
400f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
401d5b485abSSatish Balay         }
402d5b485abSSatish Balay       }
403d5b485abSSatish Balay       /* Update the headers for the current IS */
404d5b485abSSatish Balay       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
405d5b485abSSatish Balay         if ((ctr_j = ctr[j])) {
406d5b485abSSatish Balay           outdat_j        = outdat[j];
407d5b485abSSatish Balay           k               = ++outdat_j[0];
408d5b485abSSatish Balay           outdat_j[2*k]   = ctr_j;
409d5b485abSSatish Balay           outdat_j[2*k-1] = i;
410d5b485abSSatish Balay         }
411d5b485abSSatish Balay       }
412d5b485abSSatish Balay       isz[i] = isz_i;
413d5b485abSSatish Balay     }
414d5b485abSSatish Balay   }
415d5b485abSSatish Balay 
416d5b485abSSatish Balay   /*  Now  post the sends */
417*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&  s_waits1 );CHKERRQ(ierr);
418d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
419d5b485abSSatish Balay     j    = pa[i];
420ca161407SBarry Smith     ierr = MPI_Isend(outdat[j],w1[j],MPI_INT,j,tag,comm,s_waits1+i);CHKERRQ(ierr);
421d5b485abSSatish Balay   }
422d5b485abSSatish Balay 
423d5b485abSSatish Balay   /* No longer need the original indices*/
424d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
425d5b485abSSatish Balay     ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr);
426d5b485abSSatish Balay   }
427606d414cSSatish Balay   ierr = PetscFree(idx);CHKERRQ(ierr);
428d5b485abSSatish Balay 
429d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
430d5b485abSSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
431d5b485abSSatish Balay   }
432d5b485abSSatish Balay 
433d5b485abSSatish Balay   /* Do Local work*/
434df36731bSBarry Smith   ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr);
435d5b485abSSatish Balay 
436d5b485abSSatish Balay   /* Receive messages*/
437d5b485abSSatish Balay   {
438d5b485abSSatish Balay     int        index;
439d5b485abSSatish Balay 
440*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&    recv_status );CHKERRQ(ierr);
441d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
442ca161407SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&index,recv_status+i);CHKERRQ(ierr);
443d5b485abSSatish Balay     }
444d5b485abSSatish Balay 
445*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&    s_status );CHKERRQ(ierr);
446ca161407SBarry Smith     ierr     = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);
447d5b485abSSatish Balay   }
448d5b485abSSatish Balay 
449d5b485abSSatish Balay   /* Phase 1 sends are complete - deallocate buffers */
450606d414cSSatish Balay   ierr = PetscFree(outdat);CHKERRQ(ierr);
451606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
452d5b485abSSatish Balay 
453*b0a32e0cSBarry Smith   xdata = (int **)PetscMalloc((nrqr+1)*sizeof(int *));CHKERRQ(ierr);
454*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(int),&  isz1  );CHKERRQ(ierr);
455df36731bSBarry Smith   ierr  = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr);
456606d414cSSatish Balay   ierr  = PetscFree(rbuf);CHKERRQ(ierr);
457d5b485abSSatish Balay 
458d5b485abSSatish Balay   /* Send the data back*/
459d5b485abSSatish Balay   /* Do a global reduction to know the buffer space req for incoming messages*/
460d5b485abSSatish Balay   {
461d5b485abSSatish Balay     int *rw1,*rw2;
462d5b485abSSatish Balay 
463*b0a32e0cSBarry Smith ierr = PetscMalloc(2*size*sizeof(int),&(    rw1  ));CHKERRQ(ierr);
464549d3d68SSatish Balay     ierr = PetscMemzero(rw1,2*size*sizeof(int));CHKERRQ(ierr);
465d5b485abSSatish Balay     rw2  = rw1+size;
466d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
467d5b485abSSatish Balay       proc      = recv_status[i].MPI_SOURCE;
468d5b485abSSatish Balay       rw1[proc] = isz1[i];
469d5b485abSSatish Balay     }
470d5b485abSSatish Balay 
471ca161407SBarry Smith     ierr   = MPI_Allreduce(rw1,rw2,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
472d5b485abSSatish Balay     bsz1   = rw2[rank];
473606d414cSSatish Balay     ierr = PetscFree(rw1);CHKERRQ(ierr);
474d5b485abSSatish Balay   }
475d5b485abSSatish Balay 
476d5b485abSSatish Balay   /* Allocate buffers*/
477d5b485abSSatish Balay 
478d5b485abSSatish Balay   /* Allocate memory for recv buffers. Prob none if nrqr = 0 ???? */
479d5b485abSSatish Balay   len      = (nrqs+1)*sizeof(int*) + nrqs*bsz1*sizeof(int);
480*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  rbuf2    ));CHKERRQ(ierr);
481d5b485abSSatish Balay   rbuf2[0] = (int*)(rbuf2 + nrqs);
482d5b485abSSatish Balay   for (i=1; i<nrqs; ++i) rbuf2[i] = rbuf2[i-1] + bsz1;
483d5b485abSSatish Balay 
484d5b485abSSatish Balay   /* Post the receives */
485*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&  r_waits2 );CHKERRQ(ierr);
486d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
487ca161407SBarry Smith     ierr = MPI_Irecv(rbuf2[i],bsz1,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits2+i);CHKERRQ(ierr);
488d5b485abSSatish Balay   }
489d5b485abSSatish Balay 
490d5b485abSSatish Balay   /*  Now  post the sends */
491*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&  s_waits2 );CHKERRQ(ierr);
492d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
493d5b485abSSatish Balay     j    = recv_status[i].MPI_SOURCE;
494ca161407SBarry Smith     ierr = MPI_Isend(xdata[i],isz1[i],MPI_INT,j,tag,comm,s_waits2+i);CHKERRQ(ierr);
495d5b485abSSatish Balay   }
496d5b485abSSatish Balay 
497d5b485abSSatish Balay   /* receive work done on other processors*/
498d5b485abSSatish Balay   {
499d5b485abSSatish Balay     int         index,is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax;
500f1af5d2fSBarry Smith     PetscBT     table_i;
501d5b485abSSatish Balay     MPI_Status  *status2;
502d5b485abSSatish Balay 
503*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&    status2 );CHKERRQ(ierr);
504d5b485abSSatish Balay 
505d5b485abSSatish Balay     for (i=0; i<nrqs; ++i) {
506ca161407SBarry Smith       ierr = MPI_Waitany(nrqs,r_waits2,&index,status2+i);CHKERRQ(ierr);
507d5b485abSSatish Balay       /* Process the message*/
508d5b485abSSatish Balay       rbuf2_i = rbuf2[index];
509d5b485abSSatish Balay       ct1     = 2*rbuf2_i[0]+1;
510d5b485abSSatish Balay       jmax    = rbuf2[index][0];
511d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
512d5b485abSSatish Balay         max     = rbuf2_i[2*j];
513d5b485abSSatish Balay         is_no   = rbuf2_i[2*j-1];
514d5b485abSSatish Balay         isz_i   = isz[is_no];
515d5b485abSSatish Balay         data_i  = data[is_no];
516d5b485abSSatish Balay         table_i = table[is_no];
517d5b485abSSatish Balay         for (k=0; k<max; k++,ct1++) {
518d5b485abSSatish Balay           row = rbuf2_i[ct1];
519f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;}
520d5b485abSSatish Balay         }
521d5b485abSSatish Balay         isz[is_no] = isz_i;
522d5b485abSSatish Balay       }
523d5b485abSSatish Balay     }
524ca161407SBarry Smith     ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);
525606d414cSSatish Balay     ierr = PetscFree(status2);CHKERRQ(ierr);
526d5b485abSSatish Balay   }
527d5b485abSSatish Balay 
528d5b485abSSatish Balay   for (i=0; i<imax; ++i) {
529029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr);
530d5b485abSSatish Balay   }
531d5b485abSSatish Balay 
532606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
533606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
534606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
535606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
536606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
537606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
538606d414cSSatish Balay   ierr = PetscFree(table);CHKERRQ(ierr);
539606d414cSSatish Balay   ierr = PetscFree(s_status);CHKERRQ(ierr);
540606d414cSSatish Balay   ierr = PetscFree(recv_status);CHKERRQ(ierr);
541606d414cSSatish Balay   ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
542606d414cSSatish Balay   ierr = PetscFree(xdata);CHKERRQ(ierr);
543606d414cSSatish Balay   ierr = PetscFree(isz1);CHKERRQ(ierr);
5443a40ed3dSBarry Smith   PetscFunctionReturn(0);
545d5b485abSSatish Balay }
546d5b485abSSatish Balay 
5475615d1e5SSatish Balay #undef __FUNC__
548*b0a32e0cSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Local"
549d5b485abSSatish Balay /*
550df36731bSBarry Smith    MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do
551d5b485abSSatish Balay        the work on the local processor.
552d5b485abSSatish Balay 
553d5b485abSSatish Balay      Inputs:
554df36731bSBarry Smith       C      - MAT_MPIBAIJ;
555d5b485abSSatish Balay       imax - total no of index sets processed at a time;
556df36731bSBarry Smith       table  - an array of char - size = Mbs bits.
557d5b485abSSatish Balay 
558d5b485abSSatish Balay      Output:
559d5b485abSSatish Balay       isz    - array containing the count of the solution elements correspondign
560d5b485abSSatish Balay                to each index set;
561d5b485abSSatish Balay       data   - pointer to the solutions
562d5b485abSSatish Balay */
563f1af5d2fSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Local(Mat C,int imax,PetscBT *table,int *isz,int **data)
564d5b485abSSatish Balay {
565df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
566d5b485abSSatish Balay   Mat         A = c->A,B = c->B;
567df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
568df36731bSBarry Smith   int         start,end,val,max,rstart,cstart,*ai,*aj;
569d5b485abSSatish Balay   int         *bi,*bj,*garray,i,j,k,row,*data_i,isz_i;
570f1af5d2fSBarry Smith   PetscBT     table_i;
571d5b485abSSatish Balay 
5723a40ed3dSBarry Smith   PetscFunctionBegin;
573d5b485abSSatish Balay   rstart = c->rstart;
574d5b485abSSatish Balay   cstart = c->cstart;
575d5b485abSSatish Balay   ai     = a->i;
576df36731bSBarry Smith   aj     = a->j;
577d5b485abSSatish Balay   bi     = b->i;
578df36731bSBarry Smith   bj     = b->j;
579d5b485abSSatish Balay   garray = c->garray;
580d5b485abSSatish Balay 
581d5b485abSSatish Balay 
582d5b485abSSatish Balay   for (i=0; i<imax; i++) {
583d5b485abSSatish Balay     data_i  = data[i];
584d5b485abSSatish Balay     table_i = table[i];
585d5b485abSSatish Balay     isz_i   = isz[i];
586d5b485abSSatish Balay     for (j=0,max=isz[i]; j<max; j++) {
587d5b485abSSatish Balay       row   = data_i[j] - rstart;
588d5b485abSSatish Balay       start = ai[row];
589d5b485abSSatish Balay       end   = ai[row+1];
590d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Amat */
591df36731bSBarry Smith         val = aj[k] + cstart;
592f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
593d5b485abSSatish Balay       }
594d5b485abSSatish Balay       start = bi[row];
595d5b485abSSatish Balay       end   = bi[row+1];
596d5b485abSSatish Balay       for (k=start; k<end; k++) { /* Bmat */
597df36731bSBarry Smith         val = garray[bj[k]];
598f1af5d2fSBarry Smith         if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;}
599d5b485abSSatish Balay       }
600d5b485abSSatish Balay     }
601d5b485abSSatish Balay     isz[i] = isz_i;
602d5b485abSSatish Balay   }
6033a40ed3dSBarry Smith   PetscFunctionReturn(0);
604d5b485abSSatish Balay }
6055615d1e5SSatish Balay #undef __FUNC__
606*b0a32e0cSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Receive"
607d5b485abSSatish Balay /*
608df36731bSBarry Smith       MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages,
609d5b485abSSatish Balay          and return the output
610d5b485abSSatish Balay 
611d5b485abSSatish Balay          Input:
612d5b485abSSatish Balay            C    - the matrix
613d5b485abSSatish Balay            nrqr - no of messages being processed.
614d5b485abSSatish Balay            rbuf - an array of pointers to the recieved requests
615d5b485abSSatish Balay 
616d5b485abSSatish Balay          Output:
617d5b485abSSatish Balay            xdata - array of messages to be sent back
618d5b485abSSatish Balay            isz1  - size of each message
619d5b485abSSatish Balay 
620d5b485abSSatish Balay   For better efficiency perhaps we should malloc seperately each xdata[i],
621d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row
622d5b485abSSatish Balay rather then all previous rows as it is now where a single large chunck of
623d5b485abSSatish Balay memory is used.
624d5b485abSSatish Balay 
625d5b485abSSatish Balay */
626701b8814SBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,int nrqr,int **rbuf,int **xdata,int * isz1)
627d5b485abSSatish Balay {
628df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
629d5b485abSSatish Balay   Mat         A = c->A,B = c->B;
630df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
631df36731bSBarry Smith   int         rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
632d5b485abSSatish Balay   int         row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
633df36731bSBarry Smith   int         val,max1,max2,rank,Mbs,no_malloc =0,*tmp,new_estimate,ctr;
634bbd702dbSSatish Balay   int         *rbuf_i,kmax,rbuf_0,ierr;
635f1af5d2fSBarry Smith   PetscBT     xtable;
636d5b485abSSatish Balay 
6373a40ed3dSBarry Smith   PetscFunctionBegin;
638d5b485abSSatish Balay   rank   = c->rank;
639df36731bSBarry Smith   Mbs    = c->Mbs;
640d5b485abSSatish Balay   rstart = c->rstart;
641d5b485abSSatish Balay   cstart = c->cstart;
642d5b485abSSatish Balay   ai     = a->i;
643df36731bSBarry Smith   aj     = a->j;
644d5b485abSSatish Balay   bi     = b->i;
645df36731bSBarry Smith   bj     = b->j;
646d5b485abSSatish Balay   garray = c->garray;
647d5b485abSSatish Balay 
648d5b485abSSatish Balay 
649d5b485abSSatish Balay   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
650d5b485abSSatish Balay     rbuf_i  =  rbuf[i];
651d5b485abSSatish Balay     rbuf_0  =  rbuf_i[0];
652d5b485abSSatish Balay     ct     += rbuf_0;
653d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; }
654d5b485abSSatish Balay   }
655d5b485abSSatish Balay 
656701b8814SBarry Smith   if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs;
657701b8814SBarry Smith   else        max1 = 1;
658d5b485abSSatish Balay   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
659*b0a32e0cSBarry Smith   xdata[0]     = (int*)PetscMalloc(mem_estimate*sizeof(int));CHKERRQ(ierr);
660d5b485abSSatish Balay   ++no_malloc;
6616831982aSBarry Smith   ierr         = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr);
662549d3d68SSatish Balay   ierr         = PetscMemzero(isz1,nrqr*sizeof(int));CHKERRQ(ierr);
663d5b485abSSatish Balay 
664d5b485abSSatish Balay   ct3 = 0;
665d5b485abSSatish Balay   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
666d5b485abSSatish Balay     rbuf_i =  rbuf[i];
667d5b485abSSatish Balay     rbuf_0 =  rbuf_i[0];
668d5b485abSSatish Balay     ct1    =  2*rbuf_0+1;
669d5b485abSSatish Balay     ct2    =  ct1;
670d5b485abSSatish Balay     ct3    += ct1;
671d5b485abSSatish Balay     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
6726831982aSBarry Smith       ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr);
673d5b485abSSatish Balay       oct2 = ct2;
674d5b485abSSatish Balay       kmax = rbuf_i[2*j];
675d5b485abSSatish Balay       for (k=0; k<kmax; k++,ct1++) {
676d5b485abSSatish Balay         row = rbuf_i[ct1];
677f1af5d2fSBarry Smith         if (!PetscBTLookupSet(xtable,row)) {
678d5b485abSSatish Balay           if (!(ct3 < mem_estimate)) {
679d5b485abSSatish Balay             new_estimate = (int)(1.5*mem_estimate)+1;
680*b0a32e0cSBarry Smith             tmp          = (int*)PetscMalloc(new_estimate * sizeof(int));CHKERRQ(ierr);
681549d3d68SSatish Balay             ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr);
682606d414cSSatish Balay             ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
683d5b485abSSatish Balay             xdata[0]     = tmp;
684d5b485abSSatish Balay             mem_estimate = new_estimate; ++no_malloc;
685d5b485abSSatish Balay             for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
686d5b485abSSatish Balay           }
687d5b485abSSatish Balay           xdata[i][ct2++] = row;
688d5b485abSSatish Balay           ct3++;
689d5b485abSSatish Balay         }
690d5b485abSSatish Balay       }
691d5b485abSSatish Balay       for (k=oct2,max2=ct2; k<max2; k++)  {
692d5b485abSSatish Balay         row   = xdata[i][k] - rstart;
693d5b485abSSatish Balay         start = ai[row];
694d5b485abSSatish Balay         end   = ai[row+1];
695d5b485abSSatish Balay         for (l=start; l<end; l++) {
696df36731bSBarry Smith           val = aj[l] + cstart;
697f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
698d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
699d5b485abSSatish Balay               new_estimate = (int)(1.5*mem_estimate)+1;
700*b0a32e0cSBarry Smith               tmp          = (int*)PetscMalloc(new_estimate * sizeof(int));CHKERRQ(ierr);
701549d3d68SSatish Balay               ierr         = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr);
702606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
703d5b485abSSatish Balay               xdata[0]     = tmp;
704d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
705d5b485abSSatish Balay               for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
706d5b485abSSatish Balay             }
707d5b485abSSatish Balay             xdata[i][ct2++] = val;
708d5b485abSSatish Balay             ct3++;
709d5b485abSSatish Balay           }
710d5b485abSSatish Balay         }
711d5b485abSSatish Balay         start = bi[row];
712d5b485abSSatish Balay         end   = bi[row+1];
713d5b485abSSatish Balay         for (l=start; l<end; l++) {
714df36731bSBarry Smith           val = garray[bj[l]];
715f1af5d2fSBarry Smith           if (!PetscBTLookupSet(xtable,val)) {
716d5b485abSSatish Balay             if (!(ct3 < mem_estimate)) {
717d5b485abSSatish Balay               new_estimate = (int)(1.5*mem_estimate)+1;
718*b0a32e0cSBarry Smith               tmp          = (int*)PetscMalloc(new_estimate * sizeof(int));CHKERRQ(ierr);
719549d3d68SSatish Balay               ierr         = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr);
720606d414cSSatish Balay               ierr = PetscFree(xdata[0]);CHKERRQ(ierr);
721d5b485abSSatish Balay               xdata[0]     = tmp;
722d5b485abSSatish Balay               mem_estimate = new_estimate; ++no_malloc;
723d5b485abSSatish Balay               for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];}
724d5b485abSSatish Balay             }
725d5b485abSSatish Balay             xdata[i][ct2++] = val;
726d5b485abSSatish Balay             ct3++;
727d5b485abSSatish Balay           }
728d5b485abSSatish Balay         }
729d5b485abSSatish Balay       }
730d5b485abSSatish Balay       /* Update the header*/
731d5b485abSSatish Balay       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
732d5b485abSSatish Balay       xdata[i][2*j-1] = rbuf_i[2*j-1];
733d5b485abSSatish Balay     }
734d5b485abSSatish Balay     xdata[i][0] = rbuf_0;
735d5b485abSSatish Balay     xdata[i+1]  = xdata[i] + ct2;
736d5b485abSSatish Balay     isz1[i]     = ct2; /* size of each message */
737d5b485abSSatish Balay   }
7386831982aSBarry Smith   ierr = PetscBTDestroy(xtable);CHKERRQ(ierr);
739*b0a32e0cSBarry Smith   PetscLogInfo(0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %d bytes, required %d, no of mallocs = %d\n",rank,mem_estimate,ct3,no_malloc);
7403a40ed3dSBarry Smith   PetscFunctionReturn(0);
741d5b485abSSatish Balay }
742d5b485abSSatish Balay 
7437b2a1423SBarry Smith static int MatGetSubMatrices_MPIBAIJ_local(Mat,int,IS *,IS *,MatReuse,Mat *);
744a2feddc5SSatish Balay 
7455615d1e5SSatish Balay #undef __FUNC__
746*b0a32e0cSBarry Smith #define __FUNC__ "MatGetSubMatrices_MPIBAIJ"
747f1af5d2fSBarry Smith int MatGetSubMatrices_MPIBAIJ(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat **submat)
748d5b485abSSatish Balay {
74936f4e84dSSatish Balay   IS          *isrow_new,*iscol_new;
750cf4f527aSSatish Balay   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
751cf4f527aSSatish Balay   int         nmax,nstages_local,nstages,i,pos,max_no,ierr;
752a2feddc5SSatish Balay 
7533a40ed3dSBarry Smith   PetscFunctionBegin;
7543a40ed3dSBarry Smith   /* The compression and expansion should be avoided. Does'nt point
755a2feddc5SSatish Balay      out errors might change the indices hence buggey */
756a2feddc5SSatish Balay 
757*b0a32e0cSBarry Smith ierr = PetscMalloc(2*(ismax+1)*sizeof(IS),&  isrow_new );CHKERRQ(ierr);
75836f4e84dSSatish Balay   iscol_new = isrow_new + ismax;
759e757655aSSatish Balay   ierr = MatCompressIndicesSorted_MPIBAIJ(C,ismax,isrow,isrow_new);CHKERRQ(ierr);
760e757655aSSatish Balay   ierr = MatCompressIndicesSorted_MPIBAIJ(C,ismax,iscol,iscol_new);CHKERRQ(ierr);
761cf4f527aSSatish Balay 
762cf4f527aSSatish Balay   /* Allocate memory to hold all the submatrices */
763cf4f527aSSatish Balay   if (scall != MAT_REUSE_MATRIX) {
764*b0a32e0cSBarry Smith ierr = PetscMalloc((ismax+1)*sizeof(Mat),&    *submat );CHKERRQ(ierr);
765cf4f527aSSatish Balay   }
766cf4f527aSSatish Balay   /* Determine the number of stages through which submatrices are done */
767cf4f527aSSatish Balay   nmax          = 20*1000000 / (c->Nbs * sizeof(int));
768329f5518SBarry Smith   if (!nmax) nmax = 1;
769cf4f527aSSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
770cf4f527aSSatish Balay 
771cf4f527aSSatish Balay   /* Make sure every porcessor loops through the nstages */
772ca161407SBarry Smith   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr);
773cf4f527aSSatish Balay   for (i=0,pos=0; i<nstages; i++) {
774cf4f527aSSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
775cf4f527aSSatish Balay     else if (pos == ismax) max_no = 0;
776cf4f527aSSatish Balay     else                   max_no = ismax-pos;
777cf4f527aSSatish Balay     ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr);
778cf4f527aSSatish Balay     pos += max_no;
779cf4f527aSSatish Balay   }
78036f4e84dSSatish Balay 
78136f4e84dSSatish Balay   for (i=0; i<ismax; i++) {
782ca161407SBarry Smith     ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr);
783ca161407SBarry Smith     ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr);
78436f4e84dSSatish Balay   }
785606d414cSSatish Balay   ierr = PetscFree(isrow_new);CHKERRQ(ierr);
7863a40ed3dSBarry Smith   PetscFunctionReturn(0);
787a2feddc5SSatish Balay }
788a2feddc5SSatish Balay 
789233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
790233c2ff6SSatish Balay #undef __FUNC__
791233c2ff6SSatish Balay #define __FUNC__ "PetscGetProc"
792233c2ff6SSatish Balay int PetscGetProc(const int gid, const int numprocs, const int proc_gnode[], int *proc)
793233c2ff6SSatish Balay {
794233c2ff6SSatish Balay   int nGlobalNd = proc_gnode[numprocs];
795233c2ff6SSatish Balay   int fproc = (int) ((float)gid * (float)numprocs / (float)nGlobalNd + 0.5);
796233c2ff6SSatish Balay 
797233c2ff6SSatish Balay   PetscFunctionBegin;
79829bbc08cSBarry Smith   /* if(fproc < 0) SETERRQ(1,"fproc < 0");*/
799233c2ff6SSatish Balay   if (fproc > numprocs) fproc = numprocs;
800233c2ff6SSatish Balay   while (gid < proc_gnode[fproc] || gid >= proc_gnode[fproc+1]) {
801233c2ff6SSatish Balay     if (gid < proc_gnode[fproc]) fproc--;
802233c2ff6SSatish Balay     else                         fproc++;
803233c2ff6SSatish Balay   }
80429bbc08cSBarry Smith   /* if(fproc<0 || fproc>=numprocs) { SETERRQ(1,"fproc < 0 || fproc >= numprocs"); }*/
805233c2ff6SSatish Balay   *proc = fproc;
806233c2ff6SSatish Balay   PetscFunctionReturn(0);
807233c2ff6SSatish Balay }
808233c2ff6SSatish Balay #endif
809233c2ff6SSatish Balay 
810a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/
8115615d1e5SSatish Balay #undef __FUNC__
812*b0a32e0cSBarry Smith #define __FUNC__ "MatGetSubMatrices_MPIBAIJ_local"
8133eda8832SBarry Smith static int MatGetSubMatrices_MPIBAIJ_local(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat *submats)
814a2feddc5SSatish Balay {
815df36731bSBarry Smith   Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data;
816cf4f527aSSatish Balay   Mat         A = c->A;
817df36731bSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat;
818233c2ff6SSatish Balay   int         **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,start,end,size;
819233c2ff6SSatish Balay   int         **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,ierr,**rbuf1,row,proc;
82034e6ae18SSatish Balay   int         nrqs,msz,**ptr,index,*req_size,*ctr,*pa,*tmp,tcol,bsz,nrqr;
821233c2ff6SSatish Balay   int         **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2;
822233c2ff6SSatish Balay   int         **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i;
823233c2ff6SSatish Balay   int         len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
824233c2ff6SSatish Balay   int         bs=c->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB;
82536f4e84dSSatish Balay   int         cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark;
82634e6ae18SSatish Balay   int         *bmap = c->garray,ctmp,rstart=c->rstart,tag0,tag1,tag2,tag3;
827d5b485abSSatish Balay   MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
828d5b485abSSatish Balay   MPI_Request *r_waits4,*s_waits3,*s_waits4;
829d5b485abSSatish Balay   MPI_Status  *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
830d5b485abSSatish Balay   MPI_Status  *r_status3,*r_status4,*s_status4;
831d5b485abSSatish Balay   MPI_Comm    comm;
8323eda8832SBarry Smith   MatScalar   **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB;
8333eda8832SBarry Smith   MatScalar   *a_a=a->a,*b_a=b->a;
8340f5bd95cSBarry Smith   PetscTruth  flag;
8359956487aSSatish Balay   int *rw1;
83680d608c0SSatish Balay #if defined (PETSC_USE_CTABLE)
837233c2ff6SSatish Balay   int tt;
838233c2ff6SSatish Balay   PetscTable  *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1;
839233c2ff6SSatish Balay #else
840233c2ff6SSatish Balay   int         **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs;
841233c2ff6SSatish Balay #endif
842d5b485abSSatish Balay 
8433a40ed3dSBarry Smith   PetscFunctionBegin;
844d5b485abSSatish Balay   comm   = C->comm;
84534e6ae18SSatish Balay   tag0   = C->tag;
846d5b485abSSatish Balay   size   = c->size;
847d5b485abSSatish Balay   rank   = c->rank;
848d5b485abSSatish Balay 
84934e6ae18SSatish Balay   /* Get some new tags to keep the communication clean */
85034e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
85134e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
85234e6ae18SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr);
85334e6ae18SSatish Balay 
854d5b485abSSatish Balay   /* Check if the col indices are sorted */
855d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
856888f2ed8SSatish Balay     ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr);
85729bbc08cSBarry Smith     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
858d5b485abSSatish Balay   }
859d5b485abSSatish Balay 
860233c2ff6SSatish Balay   len    = (2*ismax+1)*(sizeof(int*)+ sizeof(int));
861233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE)
862233c2ff6SSatish Balay   len    += (Mbs+1)*sizeof(int);
863233c2ff6SSatish Balay #endif
864*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  irow   ));CHKERRQ(ierr);
865d5b485abSSatish Balay   icol   = irow + ismax;
866d5b485abSSatish Balay   nrow   = (int*)(icol + ismax);
867d5b485abSSatish Balay   ncol   = nrow + ismax;
868233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE)
869d5b485abSSatish Balay   rtable = ncol + ismax;
870d5b485abSSatish Balay   /* Create hash table for the mapping :row -> proc*/
871d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
872d5b485abSSatish Balay     jmax = c->rowners[i+1];
873d5b485abSSatish Balay     for (; j<jmax; j++) {
874d5b485abSSatish Balay       rtable[j] = i;
875d5b485abSSatish Balay     }
876d5b485abSSatish Balay   }
877233c2ff6SSatish Balay #endif
878233c2ff6SSatish Balay 
879233c2ff6SSatish Balay   for (i=0; i<ismax; i++) {
880233c2ff6SSatish Balay     ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
881233c2ff6SSatish Balay     ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
882233c2ff6SSatish Balay     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
883233c2ff6SSatish Balay     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
884233c2ff6SSatish Balay   }
885d5b485abSSatish Balay 
886d5b485abSSatish Balay   /* evaluate communication - mesg to who,length of mesg,and buffer space
887d5b485abSSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
888*b0a32e0cSBarry Smith ierr = PetscMalloc(size*4*sizeof(int),&(  w1     ));CHKERRQ(ierr); /* mesg size */
889d5b485abSSatish Balay   w2     = w1 + size;      /* if w2[i] marked, then a message to proc i*/
890d5b485abSSatish Balay   w3     = w2 + size;      /* no of IS that needs to be sent to proc i */
891d5b485abSSatish Balay   w4     = w3 + size;      /* temp work space used in determining w1, w2, w3 */
892549d3d68SSatish Balay   ierr   = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/
893d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
894549d3d68SSatish Balay     ierr   = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/
895d5b485abSSatish Balay     jmax   = nrow[i];
896d5b485abSSatish Balay     irow_i = irow[i];
897d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
898d5b485abSSatish Balay       row  = irow_i[j];
899233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
900233c2ff6SSatish Balay       ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr);
901233c2ff6SSatish Balay #else
902d5b485abSSatish Balay       proc = rtable[row];
903233c2ff6SSatish Balay #endif
904d5b485abSSatish Balay       w4[proc]++;
905d5b485abSSatish Balay     }
906d5b485abSSatish Balay     for (j=0; j<size; j++) {
907d5b485abSSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
908d5b485abSSatish Balay     }
909d5b485abSSatish Balay   }
910d5b485abSSatish Balay 
911d5b485abSSatish Balay   nrqs     = 0;              /* no of outgoing messages */
912e757655aSSatish Balay   msz      = 0;              /* total mesg length for all proc */
913d5b485abSSatish Balay   w1[rank] = 0;              /* no mesg sent to intself */
914d5b485abSSatish Balay   w3[rank] = 0;
915d5b485abSSatish Balay   for (i=0; i<size; i++) {
916d5b485abSSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
917d5b485abSSatish Balay   }
918*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int),&  pa );CHKERRQ(ierr); /*(proc -array)*/
919d5b485abSSatish Balay   for (i=0,j=0; i<size; i++) {
920d5b485abSSatish Balay     if (w1[i]) { pa[j] = i; j++; }
921d5b485abSSatish Balay   }
922d5b485abSSatish Balay 
923d5b485abSSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
924d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
925d5b485abSSatish Balay     j     = pa[i];
926d5b485abSSatish Balay     w1[j] += w2[j] + 2* w3[j];
927d5b485abSSatish Balay     msz   += w1[j];
928d5b485abSSatish Balay   }
929d5b485abSSatish Balay   /* Do a global reduction to determine how many messages to expect*/
930*b0a32e0cSBarry Smith ierr = PetscMalloc(2*size*sizeof(int),&(  rw1  ));CHKERRQ(ierr);
9316831982aSBarry Smith   ierr = MPI_Allreduce(w1,rw1,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr);
932d5b485abSSatish Balay   bsz  = rw1[rank];
9336831982aSBarry Smith   nrqr = rw1[size+rank];
934606d414cSSatish Balay   ierr = PetscFree(rw1);CHKERRQ(ierr);
935d5b485abSSatish Balay 
936d5b485abSSatish Balay   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
937d5b485abSSatish Balay   len      = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
938*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  rbuf1    ));CHKERRQ(ierr);
939d5b485abSSatish Balay   rbuf1[0] = (int*)(rbuf1 + nrqr);
940d5b485abSSatish Balay   for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
941d5b485abSSatish Balay 
942d5b485abSSatish Balay   /* Post the receives */
943*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&  r_waits1 );CHKERRQ(ierr);
944d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
945ca161407SBarry Smith     ierr = MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr);
946d5b485abSSatish Balay   }
947d5b485abSSatish Balay 
948d5b485abSSatish Balay   /* Allocate Memory for outgoing messages */
949d5b485abSSatish Balay   len      = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int);
950*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  sbuf1    ));CHKERRQ(ierr);
951d5b485abSSatish Balay   ptr      = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
952549d3d68SSatish Balay   ierr     = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr);
953d5b485abSSatish Balay   /* allocate memory for outgoing data + buf to receive the first reply */
954d5b485abSSatish Balay   tmp      = (int*)(ptr + size);
955d5b485abSSatish Balay   ctr      = tmp + 2*msz;
956d5b485abSSatish Balay 
957d5b485abSSatish Balay   {
958d5b485abSSatish Balay     int *iptr = tmp,ict = 0;
959d5b485abSSatish Balay     for (i=0; i<nrqs; i++) {
960d5b485abSSatish Balay       j         = pa[i];
961d5b485abSSatish Balay       iptr     += ict;
962d5b485abSSatish Balay       sbuf1[j]  = iptr;
963d5b485abSSatish Balay       ict       = w1[j];
964d5b485abSSatish Balay     }
965d5b485abSSatish Balay   }
966d5b485abSSatish Balay 
967d5b485abSSatish Balay   /* Form the outgoing messages */
968d5b485abSSatish Balay   /* Initialise the header space */
969d5b485abSSatish Balay   for (i=0; i<nrqs; i++) {
970d5b485abSSatish Balay     j           = pa[i];
971d5b485abSSatish Balay     sbuf1[j][0] = 0;
972549d3d68SSatish Balay     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr);
973d5b485abSSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
974d5b485abSSatish Balay   }
975d5b485abSSatish Balay 
976d5b485abSSatish Balay   /* Parse the isrow and copy data into outbuf */
977d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
978549d3d68SSatish Balay     ierr   = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr);
979d5b485abSSatish Balay     irow_i = irow[i];
980d5b485abSSatish Balay     jmax   = nrow[i];
981d5b485abSSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
982d5b485abSSatish Balay       row  = irow_i[j];
983233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
984233c2ff6SSatish Balay       ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr);
985233c2ff6SSatish Balay #else
986d5b485abSSatish Balay       proc = rtable[row];
987233c2ff6SSatish Balay #endif
988d5b485abSSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
989d5b485abSSatish Balay         ctr[proc]++;
990d5b485abSSatish Balay         *ptr[proc] = row;
991d5b485abSSatish Balay         ptr[proc]++;
992d5b485abSSatish Balay       }
993d5b485abSSatish Balay     }
994d5b485abSSatish Balay     /* Update the headers for the current IS */
995d5b485abSSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
996d5b485abSSatish Balay       if ((ctr_j = ctr[j])) {
997d5b485abSSatish Balay         sbuf1_j        = sbuf1[j];
998d5b485abSSatish Balay         k              = ++sbuf1_j[0];
999d5b485abSSatish Balay         sbuf1_j[2*k]   = ctr_j;
1000d5b485abSSatish Balay         sbuf1_j[2*k-1] = i;
1001d5b485abSSatish Balay       }
1002d5b485abSSatish Balay     }
1003d5b485abSSatish Balay   }
1004d5b485abSSatish Balay 
1005d5b485abSSatish Balay   /*  Now  post the sends */
1006*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&  s_waits1 );CHKERRQ(ierr);
1007d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1008d5b485abSSatish Balay     j = pa[i];
1009ca161407SBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
1010d5b485abSSatish Balay   }
1011d5b485abSSatish Balay 
1012d5b485abSSatish Balay   /* Post Recieves to capture the buffer size */
1013*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&  r_waits2 );CHKERRQ(ierr);
1014*b0a32e0cSBarry Smith   rbuf2    = (int**)PetscMalloc((nrqs+1)*sizeof(int *));CHKERRQ(ierr);
1015d5b485abSSatish Balay   rbuf2[0] = tmp + msz;
1016d5b485abSSatish Balay   for (i=1; i<nrqs; ++i) {
1017d5b485abSSatish Balay     j        = pa[i];
1018d5b485abSSatish Balay     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
1019d5b485abSSatish Balay   }
1020d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1021d5b485abSSatish Balay     j    = pa[i];
1022ca161407SBarry Smith     ierr = MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
1023d5b485abSSatish Balay   }
1024d5b485abSSatish Balay 
1025d5b485abSSatish Balay   /* Send to other procs the buf size they should allocate */
1026d5b485abSSatish Balay 
1027d5b485abSSatish Balay   /* Receive messages*/
1028*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&  s_waits2   );CHKERRQ(ierr);
1029*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&  r_status1  );CHKERRQ(ierr);
1030d5b485abSSatish Balay   len        = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*);
1031*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  sbuf2      ));CHKERRQ(ierr);
1032d5b485abSSatish Balay   req_size   = (int*)(sbuf2 + nrqr);
1033d5b485abSSatish Balay   req_source = req_size + nrqr;
1034d5b485abSSatish Balay 
1035d5b485abSSatish Balay   {
1036df36731bSBarry Smith     Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data;
1037ca161407SBarry Smith     int        *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i;
1038d5b485abSSatish Balay 
1039d5b485abSSatish Balay     for (i=0; i<nrqr; ++i) {
1040ca161407SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&index,r_status1+i);CHKERRQ(ierr);
1041d5b485abSSatish Balay       req_size[index] = 0;
1042d5b485abSSatish Balay       rbuf1_i         = rbuf1[index];
1043d5b485abSSatish Balay       start           = 2*rbuf1_i[0] + 1;
1044ca161407SBarry Smith       ierr = MPI_Get_count(r_status1+i,MPI_INT,&end);CHKERRQ(ierr);
1045*b0a32e0cSBarry Smith       sbuf2[index] = (int *)PetscMalloc(end*sizeof(int));CHKERRQ(ierr);
1046d5b485abSSatish Balay       sbuf2_i      = sbuf2[index];
1047d5b485abSSatish Balay       for (j=start; j<end; j++) {
1048d5b485abSSatish Balay         id               = rbuf1_i[j] - rstart;
1049d5b485abSSatish Balay         ncols            = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
1050d5b485abSSatish Balay         sbuf2_i[j]       = ncols;
1051d5b485abSSatish Balay         req_size[index] += ncols;
1052d5b485abSSatish Balay       }
1053d5b485abSSatish Balay       req_source[index] = r_status1[i].MPI_SOURCE;
1054d5b485abSSatish Balay       /* form the header */
1055d5b485abSSatish Balay       sbuf2_i[0]   = req_size[index];
1056d5b485abSSatish Balay       for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
1057ca161407SBarry Smith       ierr = MPI_Isend(sbuf2_i,end,MPI_INT,req_source[index],tag1,comm,s_waits2+i);CHKERRQ(ierr);
1058d5b485abSSatish Balay     }
1059d5b485abSSatish Balay   }
1060606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
1061606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
1062d5b485abSSatish Balay 
1063d5b485abSSatish Balay   /*  recv buffer sizes */
1064d5b485abSSatish Balay   /* Receive messages*/
1065d5b485abSSatish Balay 
1066*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int*),&  rbuf3     );CHKERRQ(ierr);
1067*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&  rbuf4     );CHKERRQ(ierr);
1068*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&  r_waits3  );CHKERRQ(ierr);
1069*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&  r_waits4  );CHKERRQ(ierr);
1070*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&  r_status2 );CHKERRQ(ierr);
1071d5b485abSSatish Balay 
1072d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1073ca161407SBarry Smith     ierr = MPI_Waitany(nrqs,r_waits2,&index,r_status2+i);CHKERRQ(ierr);
1074*b0a32e0cSBarry Smith     rbuf3[index] = (int *)PetscMalloc(rbuf2[index][0]*sizeof(int));CHKERRQ(ierr);
1075*b0a32e0cSBarry Smith     rbuf4[index] = (MatScalar *)PetscMalloc(rbuf2[index][0]*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1076ca161407SBarry Smith     ierr = MPI_Irecv(rbuf3[index],rbuf2[index][0],MPI_INT,
1077ca161407SBarry Smith                      r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+index);CHKERRQ(ierr);
10783eda8832SBarry Smith     ierr = MPI_Irecv(rbuf4[index],rbuf2[index][0]*bs2,MPIU_MATSCALAR,
1079ca161407SBarry Smith                      r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+index);CHKERRQ(ierr);
1080d5b485abSSatish Balay   }
1081606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
1082606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
1083d5b485abSSatish Balay 
1084d5b485abSSatish Balay   /* Wait on sends1 and sends2 */
1085*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&  s_status1 );CHKERRQ(ierr);
1086*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&  s_status2 );CHKERRQ(ierr);
1087d5b485abSSatish Balay 
1088ca161407SBarry Smith   ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);
1089ca161407SBarry Smith   ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);
1090606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
1091606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
1092606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
1093606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
1094d5b485abSSatish Balay 
1095d5b485abSSatish Balay   /* Now allocate buffers for a->j, and send them off */
1096*b0a32e0cSBarry Smith   sbuf_aj = (int **)PetscMalloc((nrqr+1)*sizeof(int *));CHKERRQ(ierr);
1097d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1098*b0a32e0cSBarry Smith   sbuf_aj[0] = (int*)PetscMalloc((j+1)*sizeof(int));CHKERRQ(ierr);
1099d5b485abSSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
1100d5b485abSSatish Balay 
1101*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&  s_waits3 );CHKERRQ(ierr);
1102d5b485abSSatish Balay   {
1103d5b485abSSatish Balay      for (i=0; i<nrqr; i++) {
1104d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
1105d5b485abSSatish Balay       sbuf_aj_i = sbuf_aj[i];
1106d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0] + 1;
1107d5b485abSSatish Balay       ct2       = 0;
1108d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
1109d5b485abSSatish Balay         kmax = rbuf1[i][2*j];
1110d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
1111e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
1112d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
1113d5b485abSSatish Balay           ncols  = nzA + nzB;
1114d5b485abSSatish Balay           cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
1115d5b485abSSatish Balay 
1116d5b485abSSatish Balay           /* load the column indices for this row into cols*/
1117d5b485abSSatish Balay           cols  = sbuf_aj_i + ct2;
1118d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
1119d5b485abSSatish Balay             if ((ctmp = bmap[cworkB[l]]) < cstart)  cols[l] = ctmp;
1120d5b485abSSatish Balay             else break;
1121d5b485abSSatish Balay           }
1122d5b485abSSatish Balay           imark = l;
1123d5b485abSSatish Balay           for (l=0; l<nzA; l++)   cols[imark+l] = cstart + cworkA[l];
1124d5b485abSSatish Balay           for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]];
1125d5b485abSSatish Balay           ct2 += ncols;
1126d5b485abSSatish Balay         }
1127d5b485abSSatish Balay       }
1128ca161407SBarry Smith       ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr);
1129d5b485abSSatish Balay     }
1130d5b485abSSatish Balay   }
1131*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&  r_status3 );CHKERRQ(ierr);
1132*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&  s_status3 );CHKERRQ(ierr);
1133d5b485abSSatish Balay 
1134d5b485abSSatish Balay   /* Allocate buffers for a->a, and send them off */
1135*b0a32e0cSBarry Smith   sbuf_aa = (MatScalar **)PetscMalloc((nrqr+1)*sizeof(MatScalar *));CHKERRQ(ierr);
1136d5b485abSSatish Balay   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1137*b0a32e0cSBarry Smith   sbuf_aa[0] = (MatScalar*)PetscMalloc((j+1)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1138a2feddc5SSatish Balay   for (i=1; i<nrqr; i++)  sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2;
1139d5b485abSSatish Balay 
1140*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&  s_waits4 );CHKERRQ(ierr);
1141d5b485abSSatish Balay   {
1142d5b485abSSatish Balay     for (i=0; i<nrqr; i++) {
1143d5b485abSSatish Balay       rbuf1_i   = rbuf1[i];
1144d5b485abSSatish Balay       sbuf_aa_i = sbuf_aa[i];
1145d5b485abSSatish Balay       ct1       = 2*rbuf1_i[0]+1;
1146d5b485abSSatish Balay       ct2       = 0;
1147d5b485abSSatish Balay       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
1148d5b485abSSatish Balay         kmax = rbuf1_i[2*j];
1149d5b485abSSatish Balay         for (k=0; k<kmax; k++,ct1++) {
1150e8e0db45SSatish Balay           row    = rbuf1_i[ct1] - rstart;
1151d5b485abSSatish Balay           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
1152d5b485abSSatish Balay           ncols  = nzA + nzB;
1153d5b485abSSatish Balay           cworkA = a_j + a_i[row];     cworkB = b_j + b_i[row];
1154a2feddc5SSatish Balay           vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2;
1155d5b485abSSatish Balay 
1156d5b485abSSatish Balay           /* load the column values for this row into vals*/
11575b83ace0SSatish Balay           vals  = sbuf_aa_i+ct2*bs2;
1158d5b485abSSatish Balay           for (l=0; l<nzB; l++) {
11593a40ed3dSBarry Smith             if ((bmap[cworkB[l]]) < cstart) {
11603eda8832SBarry Smith               ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
11613a40ed3dSBarry Smith             }
1162d5b485abSSatish Balay             else break;
1163d5b485abSSatish Balay           }
1164d5b485abSSatish Balay           imark = l;
11653a40ed3dSBarry Smith           for (l=0; l<nzA; l++) {
11663eda8832SBarry Smith             ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
11673a40ed3dSBarry Smith           }
11683a40ed3dSBarry Smith           for (l=imark; l<nzB; l++) {
11693eda8832SBarry Smith             ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
11703a40ed3dSBarry Smith           }
1171d5b485abSSatish Balay           ct2 += ncols;
1172d5b485abSSatish Balay         }
1173d5b485abSSatish Balay       }
11743eda8832SBarry Smith       ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr);
1175d5b485abSSatish Balay     }
1176d5b485abSSatish Balay   }
1177*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&  r_status4 );CHKERRQ(ierr);
1178*b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&  s_status4 );CHKERRQ(ierr);
1179606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
1180d5b485abSSatish Balay 
1181d5b485abSSatish Balay   /* Form the matrix */
1182d5b485abSSatish Balay   /* create col map */
1183d5b485abSSatish Balay   {
1184d5b485abSSatish Balay     int *icol_i;
1185233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1186233c2ff6SSatish Balay     /* Create row map*/
1187*b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&    colmaps );CHKERRQ(ierr);
1188233c2ff6SSatish Balay     for (i=0; i<ismax+1; i++) {
1189233c2ff6SSatish Balay       ierr = PetscTableCreate(((i<ismax) ? ncol[i] : ncol[i-1])+1,&colmaps[i]);CHKERRQ(ierr);
1190233c2ff6SSatish Balay     }
1191233c2ff6SSatish Balay #else
1192a2feddc5SSatish Balay     len     = (1+ismax)*sizeof(int*)+ ismax*c->Nbs*sizeof(int);
1193*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(    cmap    ));CHKERRQ(ierr);
1194d5b485abSSatish Balay     cmap[0] = (int *)(cmap + ismax);
1195549d3d68SSatish Balay     ierr    = PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(int));CHKERRQ(ierr);
1196a2feddc5SSatish Balay     for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; }
1197233c2ff6SSatish Balay #endif
1198d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1199d5b485abSSatish Balay       jmax   = ncol[i];
1200d5b485abSSatish Balay       icol_i = icol[i];
1201233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1202233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
1203233c2ff6SSatish Balay       for (j=0; j<jmax; j++) {
1204001aedefSBarry Smith 	ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr);
1205233c2ff6SSatish Balay       }
1206233c2ff6SSatish Balay #else
1207d5b485abSSatish Balay       cmap_i = cmap[i];
1208d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1209d5b485abSSatish Balay         cmap_i[icol_i[j]] = j+1;
1210d5b485abSSatish Balay       }
1211233c2ff6SSatish Balay #endif
1212d5b485abSSatish Balay     }
1213d5b485abSSatish Balay   }
1214d5b485abSSatish Balay 
1215d5b485abSSatish Balay   /* Create lens which is required for MatCreate... */
1216d5b485abSSatish Balay   for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
1217d5b485abSSatish Balay   len     = (1+ismax)*sizeof(int*)+ j*sizeof(int);
1218*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  lens    ));CHKERRQ(ierr);
1219d5b485abSSatish Balay   lens[0] = (int *)(lens + ismax);
1220549d3d68SSatish Balay   ierr    = PetscMemzero(lens[0],j*sizeof(int));CHKERRQ(ierr);
1221d5b485abSSatish Balay   for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
1222d5b485abSSatish Balay 
1223d5b485abSSatish Balay   /* Update lens from local data */
1224d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1225d5b485abSSatish Balay     jmax   = nrow[i];
1226233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1227233c2ff6SSatish Balay     lcol1_gcol1 = colmaps[i];
1228233c2ff6SSatish Balay #else
1229d5b485abSSatish Balay     cmap_i = cmap[i];
1230233c2ff6SSatish Balay #endif
1231d5b485abSSatish Balay     irow_i = irow[i];
1232d5b485abSSatish Balay     lens_i = lens[i];
1233d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1234d5b485abSSatish Balay       row  = irow_i[j];
1235233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1236233c2ff6SSatish Balay       ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr);
1237233c2ff6SSatish Balay #else
1238d5b485abSSatish Balay       proc = rtable[row];
1239233c2ff6SSatish Balay #endif
1240d5b485abSSatish Balay       if (proc == rank) {
124136f4e84dSSatish Balay         /* Get indices from matA and then from matB */
124236f4e84dSSatish Balay         row    = row - rstart;
124336f4e84dSSatish Balay         nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
124436f4e84dSSatish Balay         cworkA =  a_j + a_i[row]; cworkB = b_j + b_i[row];
1245233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1246233c2ff6SSatish Balay        for (k=0; k<nzA; k++) {
1247233c2ff6SSatish Balay 	 ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr);
1248233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1249233c2ff6SSatish Balay         }
1250233c2ff6SSatish Balay         for (k=0; k<nzB; k++) {
1251233c2ff6SSatish Balay 	  ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr);
1252233c2ff6SSatish Balay           if (tt) { lens_i[j]++; }
1253233c2ff6SSatish Balay         }
1254233c2ff6SSatish Balay #else
1255ca161407SBarry Smith         for (k=0; k<nzA; k++) {
125636f4e84dSSatish Balay           if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; }
1257ca161407SBarry Smith         }
1258ca161407SBarry Smith         for (k=0; k<nzB; k++) {
125936f4e84dSSatish Balay           if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; }
1260d5b485abSSatish Balay         }
1261233c2ff6SSatish Balay #endif
1262d5b485abSSatish Balay       }
1263a2feddc5SSatish Balay     }
1264ca161407SBarry Smith   }
1265233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1266d5b485abSSatish Balay   /* Create row map*/
1267*b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&  rowmaps );CHKERRQ(ierr);
1268233c2ff6SSatish Balay   for (i=0; i<ismax+1; i++){
1269233c2ff6SSatish Balay     ierr = PetscTableCreate((i<ismax) ? nrow[i] : nrow[i-1],&rowmaps[i]);CHKERRQ(ierr);
1270233c2ff6SSatish Balay   }
1271233c2ff6SSatish Balay #else
1272233c2ff6SSatish Balay   /* Create row map*/
1273233c2ff6SSatish Balay   len     = (1+ismax)*sizeof(int*)+ ismax*Mbs*sizeof(int);
1274*b0a32e0cSBarry Smith ierr = PetscMalloc(len,&(  rmap    ));CHKERRQ(ierr);
1275d5b485abSSatish Balay   rmap[0] = (int *)(rmap + ismax);
1276233c2ff6SSatish Balay   ierr    = PetscMemzero(rmap[0],ismax*Mbs*sizeof(int));CHKERRQ(ierr);
1277233c2ff6SSatish Balay   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;}
1278233c2ff6SSatish Balay #endif
1279d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1280d5b485abSSatish Balay     irow_i = irow[i];
1281d5b485abSSatish Balay     jmax   = nrow[i];
1282233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1283233c2ff6SSatish Balay     lrow1_grow1 = rowmaps[i];
1284233c2ff6SSatish Balay     for (j=0; j<jmax; j++) {
1285233c2ff6SSatish Balay       ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr);
1286233c2ff6SSatish Balay     }
1287233c2ff6SSatish Balay #else
1288233c2ff6SSatish Balay     rmap_i = rmap[i];
1289d5b485abSSatish Balay     for (j=0; j<jmax; j++) {
1290d5b485abSSatish Balay       rmap_i[irow_i[j]] = j;
1291d5b485abSSatish Balay     }
1292233c2ff6SSatish Balay #endif
1293d5b485abSSatish Balay   }
1294d5b485abSSatish Balay 
1295d5b485abSSatish Balay   /* Update lens from offproc data */
1296d5b485abSSatish Balay   {
1297d5b485abSSatish Balay     int *rbuf2_i,*rbuf3_i,*sbuf1_i;
1298d5b485abSSatish Balay 
1299d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1300ca161407SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits3,&i,r_status3+tmp2);CHKERRQ(ierr);
1301d5b485abSSatish Balay       index   = pa[i];
1302d5b485abSSatish Balay       sbuf1_i = sbuf1[index];
1303d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1304d5b485abSSatish Balay       ct1     = 2*jmax+1;
1305d5b485abSSatish Balay       ct2     = 0;
1306d5b485abSSatish Balay       rbuf2_i = rbuf2[i];
1307d5b485abSSatish Balay       rbuf3_i = rbuf3[i];
1308d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1309d5b485abSSatish Balay         is_no   = sbuf1_i[2*j-1];
1310d5b485abSSatish Balay         max1    = sbuf1_i[2*j];
1311d5b485abSSatish Balay         lens_i  = lens[is_no];
1312233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1313233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1314233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1315233c2ff6SSatish Balay #else
1316d5b485abSSatish Balay         cmap_i  = cmap[is_no];
1317d5b485abSSatish Balay 	rmap_i  = rmap[is_no];
1318233c2ff6SSatish Balay #endif
1319d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1320233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1321233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr);
1322233c2ff6SSatish Balay           row--;
132329bbc08cSBarry Smith           if(row < 0) { SETERRQ(1,"row not found in table"); }
1324233c2ff6SSatish Balay #else
1325d5b485abSSatish Balay           row  = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
1326233c2ff6SSatish Balay #endif
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,&tt);CHKERRQ(ierr);
1331233c2ff6SSatish Balay 	    if (tt) {
1332233c2ff6SSatish Balay               lens_i[row]++;
1333233c2ff6SSatish Balay             }
1334233c2ff6SSatish Balay #else
1335d5b485abSSatish Balay             if (cmap_i[rbuf3_i[ct2]]) {
1336d5b485abSSatish Balay               lens_i[row]++;
1337d5b485abSSatish Balay             }
1338233c2ff6SSatish Balay #endif
1339d5b485abSSatish Balay           }
1340d5b485abSSatish Balay         }
1341d5b485abSSatish Balay       }
1342d5b485abSSatish Balay     }
1343d5b485abSSatish Balay   }
1344606d414cSSatish Balay   ierr = PetscFree(r_status3);CHKERRQ(ierr);
1345606d414cSSatish Balay   ierr = PetscFree(r_waits3);CHKERRQ(ierr);
1346ca161407SBarry Smith   ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);
1347606d414cSSatish Balay   ierr = PetscFree(s_status3);CHKERRQ(ierr);
1348606d414cSSatish Balay   ierr = PetscFree(s_waits3);CHKERRQ(ierr);
1349d5b485abSSatish Balay 
1350d5b485abSSatish Balay   /* Create the submatrices */
1351d5b485abSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1352d5b485abSSatish Balay     /*
1353d5b485abSSatish Balay         Assumes new rows are same length as the old rows, hence bug!
1354d5b485abSSatish Balay     */
1355d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1356df36731bSBarry Smith       mat = (Mat_SeqBAIJ *)(submats[i]->data);
135736f4e84dSSatish Balay       if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || mat->bs != bs)) {
135829bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1359d5b485abSSatish Balay       }
13600f5bd95cSBarry Smith       ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(int),&flag);CHKERRQ(ierr);
13610f5bd95cSBarry Smith       if (flag == PETSC_FALSE) {
136229bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros");
1363d5b485abSSatish Balay       }
1364d5b485abSSatish Balay       /* Initial matrix as if empty */
1365549d3d68SSatish Balay       ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(int));CHKERRQ(ierr);
1366ce60de66SLois Curfman McInnes       submats[i]->factor = C->factor;
1367d5b485abSSatish Balay     }
1368ca161407SBarry Smith   } else {
1369d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1370ca161407SBarry Smith       ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,a->bs,nrow[i]*bs,ncol[i]*bs,0,lens[i],submats+i);CHKERRQ(ierr);
1371d5b485abSSatish Balay     }
1372d5b485abSSatish Balay   }
1373d5b485abSSatish Balay 
1374d5b485abSSatish Balay   /* Assemble the matrices */
1375d5b485abSSatish Balay   /* First assemble the local rows */
1376d5b485abSSatish Balay   {
137736f4e84dSSatish Balay     int       ilen_row,*imat_ilen,*imat_j,*imat_i;
13783eda8832SBarry Smith     MatScalar *imat_a;
1379d5b485abSSatish Balay 
1380d5b485abSSatish Balay     for (i=0; i<ismax; i++) {
1381df36731bSBarry Smith       mat       = (Mat_SeqBAIJ*)submats[i]->data;
1382d5b485abSSatish Balay       imat_ilen = mat->ilen;
1383d5b485abSSatish Balay       imat_j    = mat->j;
1384d5b485abSSatish Balay       imat_i    = mat->i;
1385d5b485abSSatish Balay       imat_a    = mat->a;
1386233c2ff6SSatish Balay 
1387233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1388233c2ff6SSatish Balay       lcol1_gcol1 = colmaps[i];
1389233c2ff6SSatish Balay       lrow1_grow1 = rowmaps[i];
1390233c2ff6SSatish Balay #else
1391d5b485abSSatish Balay       cmap_i    = cmap[i];
1392d5b485abSSatish Balay       rmap_i    = rmap[i];
1393233c2ff6SSatish Balay #endif
1394d5b485abSSatish Balay       irow_i    = irow[i];
1395d5b485abSSatish Balay       jmax      = nrow[i];
1396d5b485abSSatish Balay       for (j=0; j<jmax; j++) {
1397d5b485abSSatish Balay         row      = irow_i[j];
1398233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1399233c2ff6SSatish Balay 	ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr);
1400233c2ff6SSatish Balay #else
1401d5b485abSSatish Balay 	proc = rtable[row];
1402233c2ff6SSatish Balay #endif
1403d5b485abSSatish Balay         if (proc == rank) {
140436f4e84dSSatish Balay           row      = row - rstart;
140536f4e84dSSatish Balay           nzA      = a_i[row+1] - a_i[row];
140636f4e84dSSatish Balay           nzB      = b_i[row+1] - b_i[row];
140736f4e84dSSatish Balay           cworkA   = a_j + a_i[row];
140836f4e84dSSatish Balay           cworkB   = b_j + b_i[row];
140936f4e84dSSatish Balay           vworkA   = a_a + a_i[row]*bs2;
141036f4e84dSSatish Balay           vworkB   = b_a + b_i[row]*bs2;
1411233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1412233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr);
1413233c2ff6SSatish Balay           row--;
141429bbc08cSBarry Smith           if (row < 0) { SETERRQ(1,"row not found in table"); }
1415233c2ff6SSatish Balay #else
141636f4e84dSSatish Balay           row      = rmap_i[row + rstart];
1417233c2ff6SSatish Balay #endif
1418df36731bSBarry Smith           mat_i    = imat_i[row];
141936f4e84dSSatish Balay           mat_a    = imat_a + mat_i*bs2;
1420d5b485abSSatish Balay           mat_j    = imat_j + mat_i;
142136f4e84dSSatish Balay           ilen_row = imat_ilen[row];
142236f4e84dSSatish Balay 
142336f4e84dSSatish Balay           /* load the column indices for this row into cols*/
142436f4e84dSSatish Balay           for (l=0; l<nzB; l++) {
142536f4e84dSSatish Balay 	    if ((ctmp = bmap[cworkB[l]]) < cstart) {
1426233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1427233c2ff6SSatish Balay 	      ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr);
1428233c2ff6SSatish Balay 	      if (tcol) {
1429233c2ff6SSatish Balay #else
143036f4e84dSSatish Balay               if ((tcol = cmap_i[ctmp])) {
1431233c2ff6SSatish Balay #endif
1432df36731bSBarry Smith                 *mat_j++ = tcol - 1;
14333eda8832SBarry Smith                 ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1434549d3d68SSatish Balay                 mat_a   += bs2;
1435d5b485abSSatish Balay                 ilen_row++;
1436d5b485abSSatish Balay               }
1437ca161407SBarry Smith             } else break;
143836f4e84dSSatish Balay           }
143936f4e84dSSatish Balay           imark = l;
144036f4e84dSSatish Balay           for (l=0; l<nzA; l++) {
1441233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1442233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr);
1443233c2ff6SSatish Balay 	    if (tcol) {
1444233c2ff6SSatish Balay #else
144536f4e84dSSatish Balay 	    if ((tcol = cmap_i[cstart + cworkA[l]])) {
1446233c2ff6SSatish Balay #endif
144736f4e84dSSatish Balay               *mat_j++ = tcol - 1;
14483eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1449549d3d68SSatish Balay               mat_a   += bs2;
145036f4e84dSSatish Balay               ilen_row++;
145136f4e84dSSatish Balay             }
145236f4e84dSSatish Balay           }
145336f4e84dSSatish Balay           for (l=imark; l<nzB; l++) {
1454233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1455233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr);
1456233c2ff6SSatish Balay 	    if (tcol) {
1457233c2ff6SSatish Balay #else
145836f4e84dSSatish Balay             if ((tcol = cmap_i[bmap[cworkB[l]]])) {
1459233c2ff6SSatish Balay #endif
146036f4e84dSSatish Balay               *mat_j++ = tcol - 1;
14613eda8832SBarry Smith               ierr     = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1462549d3d68SSatish Balay               mat_a   += bs2;
146336f4e84dSSatish Balay               ilen_row++;
146436f4e84dSSatish Balay             }
146536f4e84dSSatish Balay           }
1466d5b485abSSatish Balay           imat_ilen[row] = ilen_row;
1467d5b485abSSatish Balay         }
1468d5b485abSSatish Balay       }
146936f4e84dSSatish Balay 
1470d5b485abSSatish Balay     }
1471d5b485abSSatish Balay   }
1472d5b485abSSatish Balay 
1473d5b485abSSatish Balay   /*   Now assemble the off proc rows*/
1474d5b485abSSatish Balay   {
1475d5b485abSSatish Balay     int       *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
1476d5b485abSSatish Balay     int       *imat_j,*imat_i;
14773eda8832SBarry Smith     MatScalar *imat_a,*rbuf4_i;
1478d5b485abSSatish Balay 
1479d5b485abSSatish Balay     for (tmp2=0; tmp2<nrqs; tmp2++) {
1480ca161407SBarry Smith       ierr    = MPI_Waitany(nrqs,r_waits4,&i,r_status4+tmp2);CHKERRQ(ierr);
1481d5b485abSSatish Balay       index   = pa[i];
1482d5b485abSSatish Balay       sbuf1_i = sbuf1[index];
1483d5b485abSSatish Balay       jmax    = sbuf1_i[0];
1484d5b485abSSatish Balay       ct1     = 2*jmax + 1;
1485d5b485abSSatish Balay       ct2     = 0;
1486d5b485abSSatish Balay       rbuf2_i = rbuf2[i];
1487d5b485abSSatish Balay       rbuf3_i = rbuf3[i];
1488d5b485abSSatish Balay       rbuf4_i = rbuf4[i];
1489d5b485abSSatish Balay       for (j=1; j<=jmax; j++) {
1490d5b485abSSatish Balay         is_no     = sbuf1_i[2*j-1];
1491233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1492233c2ff6SSatish Balay 	lrow1_grow1 = rowmaps[is_no];
1493233c2ff6SSatish Balay 	lcol1_gcol1 = colmaps[is_no];
1494233c2ff6SSatish Balay #else
1495d5b485abSSatish Balay         rmap_i    = rmap[is_no];
1496d5b485abSSatish Balay         cmap_i    = cmap[is_no];
1497233c2ff6SSatish Balay #endif
1498df36731bSBarry Smith         mat       = (Mat_SeqBAIJ*)submats[is_no]->data;
1499d5b485abSSatish Balay         imat_ilen = mat->ilen;
1500d5b485abSSatish Balay         imat_j    = mat->j;
1501d5b485abSSatish Balay         imat_i    = mat->i;
1502d5b485abSSatish Balay         imat_a    = mat->a;
1503d5b485abSSatish Balay         max1      = sbuf1_i[2*j];
1504d5b485abSSatish Balay         for (k=0; k<max1; k++,ct1++) {
1505d5b485abSSatish Balay           row   = sbuf1_i[ct1];
1506233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1507233c2ff6SSatish Balay 	  ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr);
1508233c2ff6SSatish Balay           row--;
150929bbc08cSBarry Smith           if(row < 0) { SETERRQ(1,"row not found in table"); }
1510233c2ff6SSatish Balay #else
1511d5b485abSSatish Balay           row   = rmap_i[row];
1512233c2ff6SSatish Balay #endif
1513d5b485abSSatish Balay           ilen  = imat_ilen[row];
1514df36731bSBarry Smith           mat_i = imat_i[row];
151536f4e84dSSatish Balay           mat_a = imat_a + mat_i*bs2;
1516d5b485abSSatish Balay           mat_j = imat_j + mat_i;
1517d5b485abSSatish Balay           max2 = rbuf2_i[ct1];
1518d5b485abSSatish Balay           for (l=0; l<max2; l++,ct2++) {
1519233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1520233c2ff6SSatish Balay 	    ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr);
1521233c2ff6SSatish Balay 	    if (tcol) {
1522233c2ff6SSatish Balay #else
1523d5b485abSSatish Balay 	    if ((tcol = cmap_i[rbuf3_i[ct2]])) {
1524233c2ff6SSatish Balay #endif
1525df36731bSBarry Smith               *mat_j++    = tcol - 1;
152636f4e84dSSatish Balay               /* *mat_a++= rbuf4_i[ct2]; */
15273eda8832SBarry Smith               ierr        = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1528549d3d68SSatish Balay               mat_a      += bs2;
1529d5b485abSSatish Balay               ilen++;
1530d5b485abSSatish Balay             }
1531d5b485abSSatish Balay           }
1532d5b485abSSatish Balay           imat_ilen[row] = ilen;
1533d5b485abSSatish Balay         }
1534d5b485abSSatish Balay       }
1535d5b485abSSatish Balay     }
1536d5b485abSSatish Balay   }
1537606d414cSSatish Balay   ierr = PetscFree(r_status4);CHKERRQ(ierr);
1538606d414cSSatish Balay   ierr = PetscFree(r_waits4);CHKERRQ(ierr);
1539ca161407SBarry Smith   ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);
1540606d414cSSatish Balay   ierr = PetscFree(s_waits4);CHKERRQ(ierr);
1541606d414cSSatish Balay   ierr = PetscFree(s_status4);CHKERRQ(ierr);
1542d5b485abSSatish Balay 
1543d5b485abSSatish Balay   /* Restore the indices */
1544d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
1545d5b485abSSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
1546d5b485abSSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
1547d5b485abSSatish Balay   }
1548d5b485abSSatish Balay 
1549d5b485abSSatish Balay   /* Destroy allocated memory */
1550606d414cSSatish Balay   ierr = PetscFree(irow);CHKERRQ(ierr);
1551606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
1552606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
1553d5b485abSSatish Balay 
1554606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
1555606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
1556d5b485abSSatish Balay   for (i=0; i<nrqr; ++i) {
1557606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
1558d5b485abSSatish Balay   }
1559d5b485abSSatish Balay   for (i=0; i<nrqs; ++i) {
1560606d414cSSatish Balay     ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr);
1561606d414cSSatish Balay     ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);
1562d5b485abSSatish Balay   }
1563d5b485abSSatish Balay 
1564606d414cSSatish Balay   ierr = PetscFree(sbuf2);CHKERRQ(ierr);
1565606d414cSSatish Balay   ierr = PetscFree(rbuf3);CHKERRQ(ierr);
1566606d414cSSatish Balay   ierr = PetscFree(rbuf4);CHKERRQ(ierr);
1567606d414cSSatish Balay   ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr);
1568606d414cSSatish Balay   ierr = PetscFree(sbuf_aj);CHKERRQ(ierr);
1569606d414cSSatish Balay   ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr);
1570606d414cSSatish Balay   ierr = PetscFree(sbuf_aa);CHKERRQ(ierr);
1571d5b485abSSatish Balay 
1572233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE)
1573233c2ff6SSatish Balay   for (i=0; i<ismax+1; i++){
1574233c2ff6SSatish Balay     ierr = PetscTableDelete(rowmaps[i]);CHKERRQ(ierr);
1575233c2ff6SSatish Balay     ierr = PetscTableDelete(colmaps[i]);CHKERRQ(ierr);
1576233c2ff6SSatish Balay   }
1577233c2ff6SSatish Balay   ierr = PetscFree(colmaps);CHKERRQ(ierr);
1578233c2ff6SSatish Balay   ierr = PetscFree(rowmaps);CHKERRQ(ierr);
1579233c2ff6SSatish Balay   /* Mark Adams */
1580233c2ff6SSatish Balay #else
1581606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
1582233c2ff6SSatish Balay   ierr = PetscFree(cmap);CHKERRQ(ierr);
1583233c2ff6SSatish Balay #endif
1584606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1585d5b485abSSatish Balay 
1586d5b485abSSatish Balay   for (i=0; i<ismax; i++) {
158736f4e84dSSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158836f4e84dSSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1589d5b485abSSatish Balay   }
159034e6ae18SSatish Balay 
15913a40ed3dSBarry Smith   PetscFunctionReturn(0);
1592d5b485abSSatish Balay }
1593ca161407SBarry Smith 
1594