xref: /petsc/src/mat/impls/sbaij/mpi/sbaijov.c (revision 5483b11dfd67de01cae3704bdb6f9c5d29f1b185)
1632d0f97SHong Zhang /*$Id: sbaijov.c,v 1.65 2001/08/06 21:15:42 bsmith Exp $*/
2632d0f97SHong Zhang 
3632d0f97SHong Zhang /*
4632d0f97SHong Zhang    Routines to compute overlapping regions of a parallel MPI matrix.
5632d0f97SHong Zhang    Used for finding submatrices that were shared across processors.
6632d0f97SHong Zhang */
7632d0f97SHong Zhang #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
8632d0f97SHong Zhang #include "petscbt.h"
9632d0f97SHong Zhang 
10632d0f97SHong Zhang static int MatIncreaseOverlap_MPISBAIJ_Once(Mat,int,IS *);
11bfc6803cSHong Zhang static int MatIncreaseOverlap_MPISBAIJ_Local(Mat,int *,int,int **,PetscBT*);
12632d0f97SHong Zhang 
13632d0f97SHong Zhang #undef __FUNCT__
14632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ"
15c910923dSHong Zhang int MatIncreaseOverlap_MPISBAIJ(Mat C,int is_max,IS is[],int ov)
16632d0f97SHong Zhang {
17d9489beaSHong Zhang   Mat_MPISBAIJ  *c = (Mat_MPISBAIJ*)C->data;
18d9489beaSHong Zhang   int           i,ierr,N=C->N, bs=c->bs;
19632d0f97SHong Zhang   IS            *is_new;
20632d0f97SHong Zhang 
21632d0f97SHong Zhang   PetscFunctionBegin;
22c910923dSHong Zhang   ierr = PetscMalloc(is_max*sizeof(IS),&is_new);CHKERRQ(ierr);
23632d0f97SHong Zhang   /* Convert the indices into block format */
2427f478b1SHong Zhang   ierr = ISCompressIndicesGeneral(N,bs,is_max,is,is_new);CHKERRQ(ierr);
25632d0f97SHong Zhang   if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
26632d0f97SHong Zhang   for (i=0; i<ov; ++i) {
27c910923dSHong Zhang     ierr = MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr);
28632d0f97SHong Zhang   }
29c910923dSHong Zhang   for (i=0; i<is_max; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);}
3027f478b1SHong Zhang   ierr = ISExpandIndicesGeneral(N,bs,is_max,is_new,is);CHKERRQ(ierr);
31c910923dSHong Zhang   for (i=0; i<is_max; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);}
32632d0f97SHong Zhang   ierr = PetscFree(is_new);CHKERRQ(ierr);
33632d0f97SHong Zhang   PetscFunctionReturn(0);
34632d0f97SHong Zhang }
35632d0f97SHong Zhang 
364a69c536SBarry Smith typedef enum {MINE,OTHER} WhoseOwner;
370472cc68SHong Zhang /*  data1, odata1 and odata2 are packed in the format (for communication):
38a2a9f22aSHong Zhang        data[0]          = is_max, no of is
39a2a9f22aSHong Zhang        data[1]          = size of is[0]
40a2a9f22aSHong Zhang         ...
41a2a9f22aSHong Zhang        data[is_max]     = size of is[is_max-1]
42a2a9f22aSHong Zhang        data[is_max + 1] = data(is[0])
43a2a9f22aSHong Zhang         ...
44a2a9f22aSHong Zhang        data[is_max+1+sum(size of is[k]), k=0,...,i-1] = data(is[i])
45a2a9f22aSHong Zhang         ...
460472cc68SHong Zhang    data2 is packed in the format (for creating output is[]):
470472cc68SHong Zhang        data[0]          = is_max, no of is
480472cc68SHong Zhang        data[1]          = size of is[0]
490472cc68SHong Zhang         ...
500472cc68SHong Zhang        data[is_max]     = size of is[is_max-1]
510472cc68SHong Zhang        data[is_max + 1] = data(is[0])
520472cc68SHong Zhang         ...
530472cc68SHong Zhang        data[is_max + 1 + Mbs*i) = data(is[i])
540472cc68SHong Zhang         ...
55a2a9f22aSHong Zhang */
56*5483b11dSHong Zhang #undef USEBcol /* not working yet! e.g., mpirun -np 3 ex92 -mat_size 100 -nd 1 (02/24/2004) */
57632d0f97SHong Zhang #undef __FUNCT__
58632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Once"
59c910923dSHong Zhang static int MatIncreaseOverlap_MPISBAIJ_Once(Mat C,int is_max,IS is[])
60632d0f97SHong Zhang {
61632d0f97SHong Zhang   Mat_MPISBAIJ  *c = (Mat_MPISBAIJ*)C->data;
62*5483b11dSHong Zhang   int         len,idx,*idx_i,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i,
630472cc68SHong Zhang               size,rank,Mbs,i,j,k,ierr,nrqs,*odata1,*odata2,
64*5483b11dSHong Zhang               tag1,tag2,flag,proc_id,**odata2_ptr,waits_idx,*ctable=0,*btable,len_b;
65*5483b11dSHong Zhang   int         nrqr,*id_r1,*len_r1,proc_end,*iwork,*len_s,*data1_i;
66bfc6803cSHong Zhang   char        *t_p;
67632d0f97SHong Zhang   MPI_Comm    comm;
68*5483b11dSHong Zhang   MPI_Request *s_waits1,*s_waits2,r_req,*r_waits1;
69632d0f97SHong Zhang   MPI_Status  *s_status,r_status;
70*5483b11dSHong Zhang   PetscBT     *table=0;  /* mark indices of this processor's is[] */
71fc70d14eSHong Zhang   PetscBT     table_i;
72bfc6803cSHong Zhang   PetscBT     otable; /* mark indices of other processors' is[] */
73*5483b11dSHong Zhang 
74632d0f97SHong Zhang   PetscFunctionBegin;
75632d0f97SHong Zhang 
76632d0f97SHong Zhang   comm = C->comm;
77632d0f97SHong Zhang   size = c->size;
78632d0f97SHong Zhang   rank = c->rank;
79632d0f97SHong Zhang   Mbs  = c->Mbs;
80632d0f97SHong Zhang 
81c910923dSHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
82c910923dSHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
83c910923dSHong Zhang 
84*5483b11dSHong Zhang   /* 1. Send this processor's is[] to other processors */
85*5483b11dSHong Zhang   /*---------------------------------------------------*/
86*5483b11dSHong Zhang   ierr = PetscMalloc(is_max*sizeof(int),&n);CHKERRQ(ierr);
87*5483b11dSHong Zhang   len = 0;
88c910923dSHong Zhang   for (i=0; i<is_max; i++) {
89632d0f97SHong Zhang     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
90632d0f97SHong Zhang     len += n[i];
91632d0f97SHong Zhang   }
92*5483b11dSHong Zhang   if (len == 0) {
93*5483b11dSHong Zhang     is_max = 0;
94*5483b11dSHong Zhang   } else {
95*5483b11dSHong Zhang     len += 1 + is_max; /* max length of data1 for one processor */
96*5483b11dSHong Zhang   }
97*5483b11dSHong Zhang   /* ierr = PetscPrintf(PETSC_COMM_SELF, " [%d] len=%d, is_max: %d\n",rank,len,is_max); */
98632d0f97SHong Zhang 
99*5483b11dSHong Zhang   ierr   = PetscMalloc(size*3*sizeof(int),&len_s);CHKERRQ(ierr);
100*5483b11dSHong Zhang   btable = len_s + size;
101*5483b11dSHong Zhang   iwork  = btable + size;
102*5483b11dSHong Zhang 
103*5483b11dSHong Zhang   ierr = PetscMalloc((size*len+1)*sizeof(int),&data1);CHKERRQ(ierr);
104*5483b11dSHong Zhang   ierr = PetscMalloc(size*sizeof(int*),&data1_start);CHKERRQ(ierr);
105*5483b11dSHong Zhang   for (i=0; i<size; i++) data1_start[i] = data1 + i*len;
106*5483b11dSHong Zhang 
107*5483b11dSHong Zhang   if (is_max){
108*5483b11dSHong Zhang     /* create hash table ctable which maps c->row to proc_id) */
109*5483b11dSHong Zhang     ierr = PetscMalloc(Mbs*sizeof(int),&ctable);CHKERRQ(ierr);
110*5483b11dSHong Zhang     for (proc_id=0,j=0; proc_id<size; proc_id++) {
111*5483b11dSHong Zhang       for (; j<c->rowners[proc_id+1]; j++) {
112*5483b11dSHong Zhang         ctable[j] = proc_id;
113*5483b11dSHong Zhang       }
114*5483b11dSHong Zhang     }
115*5483b11dSHong Zhang 
116*5483b11dSHong Zhang     /* create tables for marking indices */
117*5483b11dSHong Zhang     len_b = is_max*sizeof(PetscBT) + (Mbs/PETSC_BITS_PER_BYTE+1)*is_max*sizeof(char) + 1;
118*5483b11dSHong Zhang     ierr = PetscMalloc(len_b,&table);CHKERRQ(ierr);
119*5483b11dSHong Zhang     t_p  = (char *)(table + is_max);
120a2a9f22aSHong Zhang     for (i=0; i<is_max; i++) {
121*5483b11dSHong Zhang       table[i]  = t_p  + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
122*5483b11dSHong Zhang     }
123*5483b11dSHong Zhang 
124*5483b11dSHong Zhang #ifdef USEBcol
125*5483b11dSHong Zhang     Vec        vfrom_mpi,vto_seq;
126*5483b11dSHong Zhang     int        bs=c->bs,Bn = c->B->n,Bnbs = Bn/bs,N,low,*Bowners;
127*5483b11dSHong Zhang     PetscScalar garray[Bnbs],*array;
128*5483b11dSHong Zhang #ifdef NEW
129*5483b11dSHong Zhang     IS is_local,is_gl;
130*5483b11dSHong Zhang     ierr = ISCreateGeneral(comm,Bnbs,c->garray,&is_local);CHKERRQ(ierr);
131*5483b11dSHong Zhang     ierr = ISAllGather(is_local, &is_gl);CHKERRQ(ierr);
132*5483b11dSHong Zhang     if (rank == 1){
133*5483b11dSHong Zhang       ierr = ISView(is_gl,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
134*5483b11dSHong Zhang     }
135*5483b11dSHong Zhang     ierr = ISDestroy(is_gl);
136*5483b11dSHong Zhang     ierr = ISDestroy(is_local);
137*5483b11dSHong Zhang 
138*5483b11dSHong Zhang #endif /* NEW */
139*5483b11dSHong Zhang     /* create vfrom_mpi, a mpi vector holding local B->garray */
140*5483b11dSHong Zhang     for (i=0; i<Bnbs; i++) garray[i] = c->garray[i];
141*5483b11dSHong Zhang     ierr = VecCreateMPIWithArray(comm,Bnbs,PETSC_DECIDE,garray,&vfrom_mpi);CHKERRQ(ierr);
142*5483b11dSHong Zhang 
143*5483b11dSHong Zhang     /* create vto_seq, a seq vector for holding B->garray of all processors */
144*5483b11dSHong Zhang     ierr = VecGetSize(vfrom_mpi,&N);CHKERRQ(ierr);
145*5483b11dSHong Zhang     ierr = VecGetOwnershipRange(vfrom_mpi,&low,PETSC_NULL);CHKERRQ(ierr);
146*5483b11dSHong Zhang 
147*5483b11dSHong Zhang     ierr = PetscMalloc((size+1)*sizeof(int),&Bowners);CHKERRQ(ierr);
148*5483b11dSHong Zhang     ierr = MPI_Allgather(&low,1,MPI_INT,Bowners,1,MPI_INT,comm);CHKERRQ(ierr);
149*5483b11dSHong Zhang     Bowners[size] = N;
150*5483b11dSHong Zhang     ierr = VecConvertMPIToSeqAll(vfrom_mpi,&vto_seq);CHKERRQ(ierr);
151*5483b11dSHong Zhang 
152*5483b11dSHong Zhang     if (rank == 0){
153*5483b11dSHong Zhang       /* ierr = VecView(vto_seq,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); */
154*5483b11dSHong Zhang       for (i = 0; i<=size; i++){
155*5483b11dSHong Zhang         ierr = PetscPrintf(PETSC_COMM_SELF," Bowners[%d]=%d, N=%d\n",i,Bowners[i],N);
156*5483b11dSHong Zhang       }
157*5483b11dSHong Zhang     }
158*5483b11dSHong Zhang     ierr = VecDestroy(vfrom_mpi);CHKERRQ(ierr);
159*5483b11dSHong Zhang #endif
160*5483b11dSHong Zhang 
161*5483b11dSHong Zhang #ifdef USEBcol
162*5483b11dSHong Zhang     /* hash table table_i[idx] = 1 if idx is an index of one of the is[] array
163*5483b11dSHong Zhang                                = 0 otherwise */
164*5483b11dSHong Zhang     table_i = table[0];
165*5483b11dSHong Zhang     ierr    = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr);
166*5483b11dSHong Zhang     for (i=0; i<is_max; i++){
167a2a9f22aSHong Zhang       ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr);
168*5483b11dSHong Zhang       for (j=0; j<n[i]; j++){
169*5483b11dSHong Zhang         idx = idx_i[j];
170*5483b11dSHong Zhang         ierr = PetscBTSet(table_i,idx);CHKERRQ(ierr);
171632d0f97SHong Zhang       }
172a2a9f22aSHong Zhang       ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr);
173*5483b11dSHong Zhang     }
174*5483b11dSHong Zhang 
175*5483b11dSHong Zhang     ierr = VecGetArray(vto_seq,&array);
176*5483b11dSHong Zhang     for (i=0; i<size; i++){
177*5483b11dSHong Zhang       btable[i] = 0;
178*5483b11dSHong Zhang       for (j = Bowners[i]; j<Bowners[i+1]; j++){ /* go through B cols */
179*5483b11dSHong Zhang         idx = (int)array[j];
180*5483b11dSHong Zhang         if(PetscBTLookup(table_i,idx)){
181*5483b11dSHong Zhang           btable[i] = 1;
182*5483b11dSHong Zhang           break;
183*5483b11dSHong Zhang         }
184*5483b11dSHong Zhang       }
185*5483b11dSHong Zhang 
186*5483b11dSHong Zhang       if (!btable[i]){ /* go through A cols */
187*5483b11dSHong Zhang         for (j = c->rowners[i]; j<c->rowners[i+1]; j++){
188*5483b11dSHong Zhang           if(PetscBTLookup(table_i,j)){
189*5483b11dSHong Zhang             btable[i] = 1;
190*5483b11dSHong Zhang             break;
191*5483b11dSHong Zhang           }
192*5483b11dSHong Zhang         }
193*5483b11dSHong Zhang       }
194*5483b11dSHong Zhang     }
195*5483b11dSHong Zhang     if (rank == 1){
196*5483b11dSHong Zhang       for (i=0; i<size; i++) printf(" %d\n",btable[i]);
197*5483b11dSHong Zhang     }
198*5483b11dSHong Zhang     ierr = VecRestoreArray(vto_seq,&array);CHKERRQ(ierr);
199*5483b11dSHong Zhang     ierr = PetscFree(Bowners);CHKERRQ(ierr);
200*5483b11dSHong Zhang     ierr = VecDestroy(vto_seq);CHKERRQ(ierr);
201*5483b11dSHong Zhang #endif
202*5483b11dSHong Zhang   } /* if (is_max) */
203*5483b11dSHong Zhang 
204*5483b11dSHong Zhang   /* evaluate communication - mesg to who, length, and buffer space */
205*5483b11dSHong Zhang   for (i=0; i<size; i++) {
206*5483b11dSHong Zhang     len_s[i] = 0;
207*5483b11dSHong Zhang   }
208*5483b11dSHong Zhang 
209*5483b11dSHong Zhang   /* header of data1 */
210*5483b11dSHong Zhang   for (proc_id=0; proc_id<size; proc_id++){
211*5483b11dSHong Zhang     iwork[proc_id] = 0;
212*5483b11dSHong Zhang     *data1_start[proc_id] = is_max;
213*5483b11dSHong Zhang     data1_start[proc_id]++;
214*5483b11dSHong Zhang     for (j=0; j<is_max; j++) {
215*5483b11dSHong Zhang       if (proc_id == rank){
216*5483b11dSHong Zhang         *data1_start[proc_id] = n[j];
217*5483b11dSHong Zhang       } else {
218*5483b11dSHong Zhang         *data1_start[proc_id] = 0;
219*5483b11dSHong Zhang       }
220*5483b11dSHong Zhang       data1_start[proc_id]++;
221*5483b11dSHong Zhang     }
222*5483b11dSHong Zhang   }
223*5483b11dSHong Zhang 
224*5483b11dSHong Zhang   for (i=0; i<is_max; i++) {
225*5483b11dSHong Zhang     ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr);
226*5483b11dSHong Zhang     ierr = PetscSortInt(n[i],idx_i);CHKERRQ(ierr);
227*5483b11dSHong Zhang     for (j=0; j<n[i]; j++){
228*5483b11dSHong Zhang       idx = idx_i[j];
229*5483b11dSHong Zhang       *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */
230*5483b11dSHong Zhang       proc_end = ctable[idx];
231*5483b11dSHong Zhang       for (proc_id=0;  proc_id<=proc_end; proc_id++){ /* for others to process */
232*5483b11dSHong Zhang         if (proc_id != rank ){
233*5483b11dSHong Zhang #ifdef USEBcol
234*5483b11dSHong Zhang           if ( (proc_id < proc_end && btable[proc_id]) || proc_id == proc_end){
235*5483b11dSHong Zhang #endif
236*5483b11dSHong Zhang             *data1_start[proc_id] = idx; data1_start[proc_id]++;
237*5483b11dSHong Zhang             len_s[proc_id]++;
238*5483b11dSHong Zhang #ifdef USEBcol
239*5483b11dSHong Zhang           }
240*5483b11dSHong Zhang #endif
241*5483b11dSHong Zhang         }
242*5483b11dSHong Zhang       }
243*5483b11dSHong Zhang     }
244*5483b11dSHong Zhang 
245*5483b11dSHong Zhang     /* update header data */
246*5483b11dSHong Zhang     for (proc_id=0; proc_id<=proc_end; proc_id++){
247*5483b11dSHong Zhang       if (proc_id== rank) continue;
248*5483b11dSHong Zhang       *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id];
249*5483b11dSHong Zhang       iwork[proc_id] = len_s[proc_id] ;
250*5483b11dSHong Zhang     }
251*5483b11dSHong Zhang 
252*5483b11dSHong Zhang     ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr);
253*5483b11dSHong Zhang   } /* for (i=0; i<is_max; i++) */
254*5483b11dSHong Zhang 
255*5483b11dSHong Zhang   nrqs = 0; nrqr = 0;
256*5483b11dSHong Zhang   for (i=0; i<size; i++){
257*5483b11dSHong Zhang     data1_start[i] = data1 + i*len;
258*5483b11dSHong Zhang     if (len_s[i]){
259*5483b11dSHong Zhang       nrqs++;
260*5483b11dSHong Zhang       len_s[i] += 1 + is_max; /* add no. of header msg */
261*5483b11dSHong Zhang     }
262*5483b11dSHong Zhang   }
263*5483b11dSHong Zhang 
264*5483b11dSHong Zhang   for (i=0; i<is_max; i++) {
265a2a9f22aSHong Zhang     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
266632d0f97SHong Zhang   }
267bfc6803cSHong Zhang   ierr = PetscFree(n);CHKERRQ(ierr);
268*5483b11dSHong Zhang   if (ctable){ierr = PetscFree(ctable);CHKERRQ(ierr);}
269*5483b11dSHong Zhang 
270*5483b11dSHong Zhang   /* Determine the number of messages to expect, their lengths, from from-ids */
271*5483b11dSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&nrqr);CHKERRQ(ierr);
272*5483b11dSHong Zhang   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);CHKERRQ(ierr);
273*5483b11dSHong Zhang   /* ierr = PetscPrintf(PETSC_COMM_SELF, "[%d] nrqs: %d, nrqr: %d\n",rank,nrqs,nrqr); */
274632d0f97SHong Zhang 
275632d0f97SHong Zhang   /*  Now  post the sends */
276*5483b11dSHong Zhang   ierr = PetscMalloc(2*size*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
277*5483b11dSHong Zhang   s_waits2 = s_waits1 + size;
278632d0f97SHong Zhang   k = 0;
279*5483b11dSHong Zhang   for (proc_id=0; proc_id<size; proc_id++){  /* send data1 to processor [proc_id] */
280*5483b11dSHong Zhang     if (len_s[proc_id]){
281*5483b11dSHong Zhang       ierr = MPI_Isend(data1_start[proc_id],len_s[proc_id],MPI_INT,proc_id,tag1,comm,s_waits1+k);CHKERRQ(ierr);
282632d0f97SHong Zhang       k++;
283632d0f97SHong Zhang     }
284632d0f97SHong Zhang   }
285*5483b11dSHong Zhang   ierr = PetscFree(data1_start);CHKERRQ(ierr);
286632d0f97SHong Zhang 
287dc008846SHong Zhang   /* 2. Do local work on this processor's is[] */
288dc008846SHong Zhang   /*-------------------------------------------*/
289*5483b11dSHong Zhang   ierr = PetscBTCreate(Mbs,otable);CHKERRQ(ierr);
2900472cc68SHong Zhang   ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,data1,MINE,&data,table);CHKERRQ(ierr);
291632d0f97SHong Zhang 
292632d0f97SHong Zhang   /* 3. Receive other's is[] and process. Then send back */
293bfc6803cSHong Zhang   /*-----------------------------------------------------*/
294*5483b11dSHong Zhang   len = 0;
295*5483b11dSHong Zhang   for (i=0; i<nrqr; i++){
296*5483b11dSHong Zhang     if (len_r1[i] > len)len = len_r1[i];
297*5483b11dSHong Zhang     /* ierr = PetscPrintf(PETSC_COMM_SELF, "[%d] expect to recv len=%d from [%d]\n",rank,len_r1[i],id_r1[i]); */
298*5483b11dSHong Zhang   }
299*5483b11dSHong Zhang   ierr = PetscMalloc((len+1)*sizeof(int),&odata1);CHKERRQ(ierr);
300bfc6803cSHong Zhang   ierr = PetscMalloc(size*sizeof(int**),&odata2_ptr);CHKERRQ(ierr);
301632d0f97SHong Zhang   k = 0;
302*5483b11dSHong Zhang   while (k < nrqr){
303632d0f97SHong Zhang     /* Receive messages */
304bfc6803cSHong Zhang     ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);CHKERRQ(ierr);
305632d0f97SHong Zhang     if (flag){
306bfc6803cSHong Zhang       ierr = MPI_Get_count(&r_status,MPI_INT,&len);CHKERRQ(ierr);
307632d0f97SHong Zhang       proc_id = r_status.MPI_SOURCE;
308bfc6803cSHong Zhang       ierr = MPI_Irecv(odata1,len,MPI_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr);
309*5483b11dSHong Zhang       ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr);
310*5483b11dSHong Zhang       /* ierr = PetscPrintf(PETSC_COMM_SELF, " [%d] recv %d from [%d]\n",rank,len,proc_id); */
311632d0f97SHong Zhang 
312fc70d14eSHong Zhang       /*  Process messages */
313bfc6803cSHong Zhang       ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,&odata2_ptr[k],&otable);CHKERRQ(ierr);
314a2a9f22aSHong Zhang       odata2 = odata2_ptr[k];
315a2a9f22aSHong Zhang       len = 1 + odata2[0];
316a2a9f22aSHong Zhang       for (i=0; i<odata2[0]; i++){
317a2a9f22aSHong Zhang         len += odata2[1 + i];
318fc70d14eSHong Zhang       }
319632d0f97SHong Zhang 
320632d0f97SHong Zhang       /* Send messages back */
321*5483b11dSHong Zhang       ierr = MPI_Isend(odata2,len,MPI_INT,proc_id,tag2,comm,s_waits2+k);CHKERRQ(ierr);
322*5483b11dSHong Zhang       /* ierr = PetscPrintf(PETSC_COMM_SELF," [%d] send %d back to [%d] \n",rank,len,proc_id); */
323fc70d14eSHong Zhang       k++;
324632d0f97SHong Zhang     }
325*5483b11dSHong Zhang   }
326*5483b11dSHong Zhang   ierr = PetscFree(odata1);CHKERRQ(ierr);
327632d0f97SHong Zhang 
328fc70d14eSHong Zhang   /* 4. Receive work done on other processors, then merge */
329632d0f97SHong Zhang   /*--------------------------------------------------------*/
3300472cc68SHong Zhang   /* Allocate memory for incoming data */
331bfc6803cSHong Zhang   len = (1+is_max*(Mbs+1));
3320472cc68SHong Zhang   ierr = PetscMalloc(len*sizeof(int),&data2);CHKERRQ(ierr);
333bfc6803cSHong Zhang 
334632d0f97SHong Zhang   k = 0;
335*5483b11dSHong Zhang   while (k < nrqs){
336632d0f97SHong Zhang     /* Receive messages */
337c910923dSHong Zhang     ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status);
338632d0f97SHong Zhang     if (flag){
339bfc6803cSHong Zhang       ierr = MPI_Get_count(&r_status,MPI_INT,&len);CHKERRQ(ierr);
340632d0f97SHong Zhang       proc_id = r_status.MPI_SOURCE;
3410472cc68SHong Zhang       ierr = MPI_Irecv(data2,len,MPI_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr);
342*5483b11dSHong Zhang       ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr);
343*5483b11dSHong Zhang       /* ierr = PetscPrintf(PETSC_COMM_SELF," [%d] recv %d from [%d], data2:\n",rank,len,proc_id); */
344*5483b11dSHong Zhang       if (len > 1+is_max){ /* Add data2 into data */
3450472cc68SHong Zhang         data2_i = data2 + 1 + is_max;
346fc70d14eSHong Zhang         for (i=0; i<is_max; i++){
347fc70d14eSHong Zhang           table_i = table[i];
348bfc6803cSHong Zhang           data_i  = data + 1 + is_max + Mbs*i;
349bfc6803cSHong Zhang           isz     = data[1+i];
3500472cc68SHong Zhang           for (j=0; j<data2[1+i]; j++){
3510472cc68SHong Zhang             col = data2_i[j];
352bfc6803cSHong Zhang             if (!PetscBTLookupSet(table_i,col)) {data_i[isz++] = col;}
353fc70d14eSHong Zhang           }
354bfc6803cSHong Zhang           data[1+i] = isz;
3550472cc68SHong Zhang           if (i < is_max - 1) data2_i += data2[1+i];
356fc70d14eSHong Zhang         }
357*5483b11dSHong Zhang       }
358632d0f97SHong Zhang       k++;
359632d0f97SHong Zhang     }
360*5483b11dSHong Zhang   }
361*5483b11dSHong Zhang 
362*5483b11dSHong Zhang   /* phase 1 sends are complete */
363*5483b11dSHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&s_status);CHKERRQ(ierr);
364*5483b11dSHong Zhang   if (nrqs){
365*5483b11dSHong Zhang     ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);
366*5483b11dSHong Zhang   }
367*5483b11dSHong Zhang   ierr = PetscFree(data1);CHKERRQ(ierr);
368*5483b11dSHong Zhang 
369*5483b11dSHong Zhang   /* phase 3 sends are complete */
370*5483b11dSHong Zhang   if (nrqr){
371*5483b11dSHong Zhang     ierr = MPI_Waitall(nrqr,s_waits2,s_status);CHKERRQ(ierr);
372*5483b11dSHong Zhang   }
373*5483b11dSHong Zhang   for (k=0; k<nrqr; k++){
374*5483b11dSHong Zhang     ierr = PetscFree(odata2_ptr[k]);CHKERRQ(ierr);
375*5483b11dSHong Zhang   }
376*5483b11dSHong Zhang   ierr = PetscFree(odata2_ptr);CHKERRQ(ierr);
377632d0f97SHong Zhang 
378c910923dSHong Zhang   /* 5. Create new is[] */
379c910923dSHong Zhang   /*--------------------*/
380c910923dSHong Zhang   for (i=0; i<is_max; i++) {
381bfc6803cSHong Zhang     data_i = data + 1 + is_max + Mbs*i;
382bfc6803cSHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,is+i);CHKERRQ(ierr);
383632d0f97SHong Zhang   }
384*5483b11dSHong Zhang 
385a2a9f22aSHong Zhang   ierr = PetscFree(data2);CHKERRQ(ierr);
3860472cc68SHong Zhang   ierr = PetscFree(data);CHKERRQ(ierr);
387*5483b11dSHong Zhang   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
388632d0f97SHong Zhang   ierr = PetscFree(s_status);CHKERRQ(ierr);
389*5483b11dSHong Zhang   if (table) {ierr = PetscFree(table);CHKERRQ(ierr);}
390bfc6803cSHong Zhang   ierr = PetscBTDestroy(otable);CHKERRQ(ierr);
391*5483b11dSHong Zhang 
392*5483b11dSHong Zhang   ierr = PetscFree(len_s);CHKERRQ(ierr);
393*5483b11dSHong Zhang   ierr = PetscFree(id_r1);CHKERRQ(ierr);
394*5483b11dSHong Zhang   ierr = PetscFree(len_r1);CHKERRQ(ierr);
395*5483b11dSHong Zhang 
396632d0f97SHong Zhang   PetscFunctionReturn(0);
397632d0f97SHong Zhang }
398632d0f97SHong Zhang 
399632d0f97SHong Zhang #undef __FUNCT__
400632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Local"
401632d0f97SHong Zhang /*
402dc008846SHong Zhang    MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do
403632d0f97SHong Zhang        the work on the local processor.
404632d0f97SHong Zhang 
405632d0f97SHong Zhang      Inputs:
406632d0f97SHong Zhang       C      - MAT_MPISBAIJ;
407bfc6803cSHong Zhang       data   - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format.
408bfc6803cSHong Zhang       whose  - whose is[] to be processed,
409bfc6803cSHong Zhang                MINE:  this processor's is[]
410bfc6803cSHong Zhang                OTHER: other processor's is[]
411632d0f97SHong Zhang      Output:
4120472cc68SHong Zhang        data_new  - whose = MINE:
4130472cc68SHong Zhang                      holds input and newly found indices in the same format as data
4140472cc68SHong Zhang                    whose = OTHER:
4150472cc68SHong Zhang                      only holds the newly found indices
4160472cc68SHong Zhang        table     - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'.
417632d0f97SHong Zhang */
418bfc6803cSHong Zhang static int MatIncreaseOverlap_MPISBAIJ_Local(Mat C,int *data,int whose,int **data_new,PetscBT *table)
419632d0f97SHong Zhang {
420632d0f97SHong Zhang   Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data;
421dc008846SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(c->A)->data;
422dc008846SHong Zhang   Mat_SeqBAIJ  *b = (Mat_SeqBAIJ*)(c->B)->data;
42331f99336SHong Zhang   int          ierr,row,mbs,Mbs,*nidx,*nidx_i,col,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l;
424dc008846SHong Zhang   int          a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n;
425bfc6803cSHong Zhang   PetscBT      table0;  /* mark the indices of input is[] for look up */
426bfc6803cSHong Zhang   PetscBT      table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */
427632d0f97SHong Zhang 
428632d0f97SHong Zhang   PetscFunctionBegin;
42931f99336SHong Zhang   Mbs = c->Mbs; mbs = a->mbs;
430dc008846SHong Zhang   ai = a->i; aj = a->j;
431dc008846SHong Zhang   bi = b->i; bj = b->j;
432632d0f97SHong Zhang   garray = c->garray;
433dc008846SHong Zhang   rstart = c->rstart;
434dc008846SHong Zhang   is_max = data[0];
435632d0f97SHong Zhang 
436fc70d14eSHong Zhang   ierr = PetscBTCreate(Mbs,table0);CHKERRQ(ierr);
437fc70d14eSHong Zhang 
438fc70d14eSHong Zhang   ierr = PetscMalloc((1+is_max*(Mbs+1))*sizeof(int),&nidx);CHKERRQ(ierr);
439fc70d14eSHong Zhang   nidx[0] = is_max;
440fc70d14eSHong Zhang 
441fc70d14eSHong Zhang   idx_i  = data + is_max + 1; /* ptr to input is[0] array */
442bfc6803cSHong Zhang   nidx_i = nidx + is_max + 1; /* ptr to output is[0] array */
443dc008846SHong Zhang   for (i=0; i<is_max; i++) { /* for each is */
444dc008846SHong Zhang     isz  = 0;
445fc70d14eSHong Zhang     n = data[1+i]; /* size of input is[i] */
446dc008846SHong Zhang 
447bfc6803cSHong Zhang     /* initialize table_i, set table0 */
448bfc6803cSHong Zhang     if (whose == MINE){ /* process this processor's is[] */
449bfc6803cSHong Zhang       table_i = table[i];
4500472cc68SHong Zhang       nidx_i  = nidx + 1+ is_max + Mbs*i;
451bfc6803cSHong Zhang     } else {            /* process other processor's is[] - only use one temp table */
452bfc6803cSHong Zhang       table_i = *table;
453bfc6803cSHong Zhang     }
454bfc6803cSHong Zhang     ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr);
455bfc6803cSHong Zhang     ierr = PetscBTMemzero(Mbs,table0);CHKERRQ(ierr);
456fc70d14eSHong Zhang     if (n > 0) {
457bfc6803cSHong Zhang       isz0 = 0;
458dc008846SHong Zhang       for (j=0; j<n; j++){
459dc008846SHong Zhang         col = idx_i[j];
460dc008846SHong Zhang         if (col >= Mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"index col %d >= Mbs %d",col,Mbs);
461bfc6803cSHong Zhang         if(!PetscBTLookupSet(table_i,col)) {
462bfc6803cSHong Zhang           ierr = PetscBTSet(table0,col);CHKERRQ(ierr);
4630472cc68SHong Zhang           if (whose == MINE) {nidx_i[isz0] = col;}
464bfc6803cSHong Zhang           isz0++;
465bfc6803cSHong Zhang         }
466632d0f97SHong Zhang       }
467dc008846SHong Zhang 
4680472cc68SHong Zhang       if (whose == MINE) {isz = isz0;}
469fc70d14eSHong Zhang       k = 0;  /* no. of indices from input is[i] that have been examined */
470dc008846SHong Zhang       for (row=0; row<mbs; row++){
471dc008846SHong Zhang         a_start = ai[row]; a_end = ai[row+1];
472dc008846SHong Zhang         b_start = bi[row]; b_end = bi[row+1];
4730472cc68SHong Zhang         if (PetscBTLookup(table0,row+rstart)){ /* row is on input is[i]:
4740472cc68SHong Zhang            do row search: collect all col in this row */
475dc008846SHong Zhang           for (l = a_start; l<a_end ; l++){ /* Amat */
476dc008846SHong Zhang             col = aj[l] + rstart;
477fc70d14eSHong Zhang             if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;}
478dc008846SHong Zhang           }
479dc008846SHong Zhang           for (l = b_start; l<b_end ; l++){ /* Bmat */
480dc008846SHong Zhang             col = garray[bj[l]];
481fc70d14eSHong Zhang             if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;}
482dc008846SHong Zhang           }
483dc008846SHong Zhang           k++;
484dc008846SHong Zhang           if (k >= isz0) break; /* for (row=0; row<mbs; row++) */
4850472cc68SHong Zhang         } else { /* row is not on input is[i]:
4860472cc68SHong Zhang           do col serach: add row onto nidx_i if there is a col in nidx_i */
487dc008846SHong Zhang           for (l = a_start; l<a_end ; l++){ /* Amat */
488dc008846SHong Zhang             col = aj[l] + rstart;
489dc008846SHong Zhang             if (PetscBTLookup(table0,col)){
490fc70d14eSHong Zhang               if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;}
491dc008846SHong Zhang               break; /* for l = start; l<end ; l++) */
492632d0f97SHong Zhang             }
493632d0f97SHong Zhang           }
494dc008846SHong Zhang           for (l = b_start; l<b_end ; l++){ /* Bmat */
495dc008846SHong Zhang             col = garray[bj[l]];
496dc008846SHong Zhang             if (PetscBTLookup(table0,col)){
497fc70d14eSHong Zhang               if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;}
498dc008846SHong Zhang               break; /* for l = start; l<end ; l++) */
499632d0f97SHong Zhang             }
500dc008846SHong Zhang           }
501dc008846SHong Zhang         }
502bfc6803cSHong Zhang       }
503fc70d14eSHong Zhang     } /* if (n > 0) */
504dc008846SHong Zhang 
505dc008846SHong Zhang     if (i < is_max - 1){
506fc70d14eSHong Zhang       idx_i  += n;   /* ptr to input is[i+1] array */
507bfc6803cSHong Zhang       nidx_i += isz; /* ptr to output is[i+1] array */
508dc008846SHong Zhang     }
509fc70d14eSHong Zhang     nidx[1+i] = isz; /* size of new is[i] */
5101ab97a2aSSatish Balay   } /* for each is */
511dc008846SHong Zhang   *data_new = nidx;
512dc008846SHong Zhang   ierr = PetscBTDestroy(table0);CHKERRQ(ierr);
513632d0f97SHong Zhang 
514632d0f97SHong Zhang   PetscFunctionReturn(0);
515632d0f97SHong Zhang }
516632d0f97SHong Zhang 
517632d0f97SHong Zhang 
518