xref: /petsc/src/mat/impls/sbaij/mpi/sbaijov.c (revision be2f4de2ef64ba82ebcbfccb98c2538d8cd32109)
1 
2 /*
3    Routines to compute overlapping regions of a parallel MPI matrix.
4    Used for finding submatrices that were shared across processors.
5 */
6 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
7 #include <petscbt.h>
8 
9 static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat,PetscInt,IS*);
10 static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat,PetscInt*,PetscInt,PetscInt*,PetscBT*);
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ"
14 PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat C,PetscInt is_max,IS is[],PetscInt ov)
15 {
16   PetscErrorCode ierr;
17   PetscInt       i,N=C->cmap->N, bs=C->rmap->bs,M=C->rmap->N,Mbs=M/bs,*nidx,isz,iov;
18   IS             *is_new,*is_row;
19   Mat            *submats;
20   Mat_MPISBAIJ   *c=(Mat_MPISBAIJ*)C->data;
21   Mat_SeqSBAIJ   *asub_i;
22   PetscBT        table;
23   PetscInt       *ai,brow,nz,nis,l,nmax,nstages_local,nstages,max_no,pos;
24   const PetscInt *idx;
25 
26   PetscFunctionBegin;
27   ierr = PetscMalloc(is_max*sizeof(IS),&is_new);CHKERRQ(ierr);
28   /* Convert the indices into block format */
29   ierr = ISCompressIndicesGeneral(N,bs,is_max,is,is_new);CHKERRQ(ierr);
30   if (ov < 0){ SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
31 
32   /* ----- previous non-scalable implementation ----- */
33   PetscBool flg=PETSC_FALSE;
34   ierr = PetscOptionsHasName(PETSC_NULL, "-IncreaseOverlap_old", &flg);CHKERRQ(ierr);
35   if (flg){ /* previous non-scalable implementation */
36     printf("use previous non-scalable implementation...\n");
37     for (i=0; i<ov; ++i) {
38       ierr = MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr);
39     }
40   } else { /* scalable implementation using modified BAIJ routines */
41 
42   ierr = PetscMalloc((Mbs+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
43   ierr = PetscBTCreate(Mbs,table);CHKERRQ(ierr); /* for column search */
44 
45   /* Create is_row */
46   ierr = PetscMalloc(is_max*sizeof(IS **),&is_row);CHKERRQ(ierr);
47   ierr = ISCreateStride(PETSC_COMM_SELF,Mbs,0,1,&is_row[0]);CHKERRQ(ierr);
48   for (i=1; i<is_max; i++) is_row[i] = is_row[0]; /* reuse is_row[0] */
49 
50   /* Allocate memory to hold all the submatrices - Modified from MatGetSubMatrices_MPIBAIJ() */
51   ierr = PetscMalloc((is_max+1)*sizeof(Mat),&submats);CHKERRQ(ierr);
52 
53   /* Determine the number of stages through which submatrices are done */
54   nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt));
55   if (!nmax) nmax = 1;
56   nstages_local = is_max/nmax + ((is_max % nmax)?1:0);
57 
58   /* Make sure every processor loops through the nstages */
59   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);CHKERRQ(ierr);
60 
61   for (iov=0; iov<ov; ++iov) {
62     /* 1) Get submats for column search */
63     for (i=0,pos=0; i<nstages; i++) {
64       if (pos+nmax <= is_max) max_no = nmax;
65       else if (pos == is_max) max_no = 0;
66       else                   max_no = is_max-pos;
67 
68       c->ijonly = PETSC_TRUE;
69       ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,is_row+pos,is_new+pos,MAT_INITIAL_MATRIX,submats+pos);CHKERRQ(ierr);
70       pos += max_no;
71     }
72 
73     /* 2) Row search */
74     ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr);
75 
76     /* 3) Column search */
77     for (i=0; i<is_max; i++){
78       asub_i = (Mat_SeqSBAIJ*)submats[i]->data;
79       ai=asub_i->i;;
80 
81       /* put is_new obtained from MatIncreaseOverlap_MPIBAIJ() to table */
82       ierr = PetscBTMemzero(Mbs,table);CHKERRQ(ierr);
83 
84       ierr = ISGetIndices(is_new[i],&idx);CHKERRQ(ierr);
85       ierr = ISGetLocalSize(is_new[i],&nis);CHKERRQ(ierr);
86       for (l=0; l<nis; l++) {
87         ierr = PetscBTSet(table,idx[l]);CHKERRQ(ierr);
88         nidx[l] = idx[l];
89       }
90       isz = nis;
91 
92       /* add column entries to table */
93       for (brow=0; brow<Mbs; brow++){
94         nz = ai[brow+1] - ai[brow];
95         if (nz) {
96           if (!PetscBTLookupSet(table,brow)) nidx[isz++] = brow;
97         }
98       }
99       ierr = ISRestoreIndices(is_new[i],&idx);CHKERRQ(ierr);
100       ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);
101 
102       /* create updated is_new */
103       ierr = ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,is_new+i);CHKERRQ(ierr);
104     }
105 
106     /* Free tmp spaces */
107     ierr = PetscBTDestroy(table);CHKERRQ(ierr);
108     for (i=0; i<is_max; i++){
109       ierr = MatDestroy(&submats[i]);CHKERRQ(ierr);
110     }
111     ierr = PetscFree(submats);CHKERRQ(ierr);
112   }
113 
114   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
115   ierr = PetscFree(submats);CHKERRQ(ierr);
116   ierr = ISDestroy(&is_row[0]);CHKERRQ(ierr);
117   ierr = PetscFree(is_row);CHKERRQ(ierr);
118   ierr = PetscFree(nidx);CHKERRQ(ierr);
119 
120   }
121   //--------end of new----------
122 
123   for (i=0; i<is_max; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);}
124   ierr = ISExpandIndicesGeneral(N,bs,is_max,is_new,is);CHKERRQ(ierr);
125 
126   for (i=0; i<is_max; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);}
127   ierr = PetscFree(is_new);CHKERRQ(ierr);
128   PetscFunctionReturn(0);
129 }
130 
131 typedef enum {MINE,OTHER} WhoseOwner;
132 /*  data1, odata1 and odata2 are packed in the format (for communication):
133        data[0]          = is_max, no of is
134        data[1]          = size of is[0]
135         ...
136        data[is_max]     = size of is[is_max-1]
137        data[is_max + 1] = data(is[0])
138         ...
139        data[is_max+1+sum(size of is[k]), k=0,...,i-1] = data(is[i])
140         ...
141    data2 is packed in the format (for creating output is[]):
142        data[0]          = is_max, no of is
143        data[1]          = size of is[0]
144         ...
145        data[is_max]     = size of is[is_max-1]
146        data[is_max + 1] = data(is[0])
147         ...
148        data[is_max + 1 + Mbs*i) = data(is[i])
149         ...
150 */
151 #undef __FUNCT__
152 #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Once"
153 static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat C,PetscInt is_max,IS is[])
154 {
155   Mat_MPISBAIJ  *c = (Mat_MPISBAIJ*)C->data;
156   PetscErrorCode ierr;
157   PetscMPIInt    size,rank,tag1,tag2,*len_s,nrqr,nrqs,*id_r1,*len_r1,flag,len;
158   const PetscInt *idx_i;
159   PetscInt       idx,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i,
160                  Mbs,i,j,k,*odata1,*odata2,
161                  proc_id,**odata2_ptr,*ctable=0,*btable,len_max,len_est;
162   PetscInt       proc_end=0,*iwork,len_unused,nodata2;
163   PetscInt       ois_max; /* max no of is[] in each of processor */
164   char           *t_p;
165   MPI_Comm       comm;
166   MPI_Request    *s_waits1,*s_waits2,r_req;
167   MPI_Status     *s_status,r_status;
168   PetscBT        *table;  /* mark indices of this processor's is[] */
169   PetscBT        table_i;
170   PetscBT        otable; /* mark indices of other processors' is[] */
171   PetscInt       bs=C->rmap->bs,Bn = c->B->cmap->n,Bnbs = Bn/bs,*Bowners;
172   IS             garray_local,garray_gl;
173 
174   PetscFunctionBegin;
175   comm = ((PetscObject)C)->comm;
176   size = c->size;
177   rank = c->rank;
178   Mbs  = c->Mbs;
179 
180   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
181   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
182 
183   /* create tables used in
184      step 1: table[i] - mark c->garray of proc [i]
185      step 3: table[i] - mark indices of is[i] when whose=MINE
186              table[0] - mark incideces of is[] when whose=OTHER */
187   len = PetscMax(is_max, size);CHKERRQ(ierr);
188   ierr = PetscMalloc2(len,PetscBT,&table,(Mbs/PETSC_BITS_PER_BYTE+1)*len,char,&t_p);CHKERRQ(ierr);
189   for (i=0; i<len; i++) {
190     table[i]  = t_p  + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
191   }
192 
193   ierr = MPI_Allreduce(&is_max,&ois_max,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
194 
195   /* 1. Send this processor's is[] to other processors */
196   /*---------------------------------------------------*/
197   /* allocate spaces */
198   ierr = PetscMalloc(is_max*sizeof(PetscInt),&n);CHKERRQ(ierr);
199   len = 0;
200   for (i=0; i<is_max; i++) {
201     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
202     len += n[i];
203   }
204   if (!len) {
205     is_max = 0;
206   } else {
207     len += 1 + is_max; /* max length of data1 for one processor */
208   }
209 
210 
211   ierr = PetscMalloc((size*len+1)*sizeof(PetscInt),&data1);CHKERRQ(ierr);
212   ierr = PetscMalloc(size*sizeof(PetscInt*),&data1_start);CHKERRQ(ierr);
213   for (i=0; i<size; i++) data1_start[i] = data1 + i*len;
214 
215   ierr = PetscMalloc4(size,PetscInt,&len_s,size,PetscInt,&btable,size,PetscInt,&iwork,size+1,PetscInt,&Bowners);CHKERRQ(ierr);
216 
217   /* gather c->garray from all processors */
218   ierr = ISCreateGeneral(comm,Bnbs,c->garray,PETSC_COPY_VALUES,&garray_local);CHKERRQ(ierr);
219   ierr = ISAllGather(garray_local, &garray_gl);CHKERRQ(ierr);
220   ierr = ISDestroy(&garray_local);CHKERRQ(ierr);
221   ierr = MPI_Allgather(&Bnbs,1,MPIU_INT,Bowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
222   Bowners[0] = 0;
223   for (i=0; i<size; i++) Bowners[i+1] += Bowners[i];
224 
225   if (is_max){
226     /* hash table ctable which maps c->row to proc_id) */
227     ierr = PetscMalloc(Mbs*sizeof(PetscInt),&ctable);CHKERRQ(ierr);
228     for (proc_id=0,j=0; proc_id<size; proc_id++) {
229       for (; j<C->rmap->range[proc_id+1]/bs; j++) {
230         ctable[j] = proc_id;
231       }
232     }
233 
234     /* hash tables marking c->garray */
235     ierr = ISGetIndices(garray_gl,&idx_i);
236     for (i=0; i<size; i++){
237       table_i = table[i];
238       ierr    = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr);
239       for (j = Bowners[i]; j<Bowners[i+1]; j++){ /* go through B cols of proc[i]*/
240         ierr = PetscBTSet(table_i,idx_i[j]);CHKERRQ(ierr);
241       }
242     }
243     ierr = ISRestoreIndices(garray_gl,&idx_i);CHKERRQ(ierr);
244   }  /* if (is_max) */
245   ierr = ISDestroy(&garray_gl);CHKERRQ(ierr);
246 
247   /* evaluate communication - mesg to who, length, and buffer space */
248   for (i=0; i<size; i++) len_s[i] = 0;
249 
250   /* header of data1 */
251   for (proc_id=0; proc_id<size; proc_id++){
252     iwork[proc_id] = 0;
253     *data1_start[proc_id] = is_max;
254     data1_start[proc_id]++;
255     for (j=0; j<is_max; j++) {
256       if (proc_id == rank){
257         *data1_start[proc_id] = n[j];
258       } else {
259         *data1_start[proc_id] = 0;
260       }
261       data1_start[proc_id]++;
262     }
263   }
264 
265   for (i=0; i<is_max; i++) {
266     ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr);
267     for (j=0; j<n[i]; j++){
268       idx = idx_i[j];
269       *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */
270       proc_end = ctable[idx];
271       for (proc_id=0;  proc_id<=proc_end; proc_id++){ /* for others to process */
272         if (proc_id == rank ) continue; /* done before this loop */
273         if (proc_id < proc_end && !PetscBTLookup(table[proc_id],idx))
274           continue;   /* no need for sending idx to [proc_id] */
275         *data1_start[proc_id] = idx; data1_start[proc_id]++;
276         len_s[proc_id]++;
277       }
278     }
279     /* update header data */
280     for (proc_id=0; proc_id<size; proc_id++){
281       if (proc_id== rank) continue;
282       *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id];
283       iwork[proc_id] = len_s[proc_id] ;
284     }
285     ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr);
286   }
287 
288   nrqs = 0; nrqr = 0;
289   for (i=0; i<size; i++){
290     data1_start[i] = data1 + i*len;
291     if (len_s[i]){
292       nrqs++;
293       len_s[i] += 1 + is_max; /* add no. of header msg */
294     }
295   }
296 
297   for (i=0; i<is_max; i++) {
298     ierr = ISDestroy(&is[i]);CHKERRQ(ierr);
299   }
300   ierr = PetscFree(n);CHKERRQ(ierr);
301   ierr = PetscFree(ctable);CHKERRQ(ierr);
302 
303   /* Determine the number of messages to expect, their lengths, from from-ids */
304   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&nrqr);CHKERRQ(ierr);
305   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);CHKERRQ(ierr);
306 
307   /*  Now  post the sends */
308   ierr = PetscMalloc2(size,MPI_Request,&s_waits1,size,MPI_Request,&s_waits2);CHKERRQ(ierr);
309   k = 0;
310   for (proc_id=0; proc_id<size; proc_id++){  /* send data1 to processor [proc_id] */
311     if (len_s[proc_id]){
312       ierr = MPI_Isend(data1_start[proc_id],len_s[proc_id],MPIU_INT,proc_id,tag1,comm,s_waits1+k);CHKERRQ(ierr);
313       k++;
314     }
315   }
316 
317   /* 2. Receive other's is[] and process. Then send back */
318   /*-----------------------------------------------------*/
319   len = 0;
320   for (i=0; i<nrqr; i++){
321     if (len_r1[i] > len)len = len_r1[i];
322   }
323   ierr = PetscFree(len_r1);CHKERRQ(ierr);
324   ierr = PetscFree(id_r1);CHKERRQ(ierr);
325 
326   for (proc_id=0; proc_id<size; proc_id++)
327     len_s[proc_id] = iwork[proc_id] = 0;
328 
329   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&odata1);CHKERRQ(ierr);
330   ierr = PetscMalloc(size*sizeof(PetscInt**),&odata2_ptr);CHKERRQ(ierr);
331   ierr = PetscBTCreate(Mbs,otable);CHKERRQ(ierr);
332 
333   len_max = ois_max*(Mbs+1);  /* max space storing all is[] for each receive */
334   len_est = 2*len_max; /* estimated space of storing is[] for all receiving messages */
335   ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr);
336   nodata2 = 0;       /* nodata2+1: num of PetscMalloc(,&odata2_ptr[]) called */
337   odata2_ptr[nodata2] = odata2;
338   len_unused = len_est; /* unused space in the array odata2_ptr[nodata2]-- needs to be >= len_max  */
339 
340   k = 0;
341   while (k < nrqr){
342     /* Receive messages */
343     ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);CHKERRQ(ierr);
344     if (flag){
345       ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRQ(ierr);
346       proc_id = r_status.MPI_SOURCE;
347       ierr = MPI_Irecv(odata1,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr);
348       ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr);
349 
350       /*  Process messages */
351       /*  make sure there is enough unused space in odata2 array */
352       if (len_unused < len_max){ /* allocate more space for odata2 */
353         ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr);
354         odata2_ptr[++nodata2] = odata2;
355         len_unused = len_est;
356       }
357 
358       ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,odata2,&otable);CHKERRQ(ierr);
359       len = 1 + odata2[0];
360       for (i=0; i<odata2[0]; i++){
361         len += odata2[1 + i];
362       }
363 
364       /* Send messages back */
365       ierr = MPI_Isend(odata2,len,MPIU_INT,proc_id,tag2,comm,s_waits2+k);CHKERRQ(ierr);
366       k++;
367       odata2     += len;
368       len_unused -= len;
369       len_s[proc_id] = len; /* num of messages sending back to [proc_id] by this proc */
370     }
371   }
372   ierr = PetscFree(odata1);CHKERRQ(ierr);
373   ierr = PetscBTDestroy(otable);CHKERRQ(ierr);
374 
375   /* 3. Do local work on this processor's is[] */
376   /*-------------------------------------------*/
377   /* make sure there is enough unused space in odata2(=data) array */
378   len_max = is_max*(Mbs+1); /* max space storing all is[] for this processor */
379   if (len_unused < len_max){ /* allocate more space for odata2 */
380     ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr);
381     odata2_ptr[++nodata2] = odata2;
382     len_unused = len_est;
383   }
384 
385   data = odata2;
386   ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,data,table);CHKERRQ(ierr);
387   ierr = PetscFree(data1_start);CHKERRQ(ierr);
388 
389   /* 4. Receive work done on other processors, then merge */
390   /*------------------------------------------------------*/
391   /* get max number of messages that this processor expects to recv */
392   ierr = MPI_Allreduce(len_s,iwork,size,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
393   ierr = PetscMalloc((iwork[rank]+1)*sizeof(PetscInt),&data2);CHKERRQ(ierr);
394   ierr = PetscFree4(len_s,btable,iwork,Bowners);CHKERRQ(ierr);
395 
396   k = 0;
397   while (k < nrqs){
398     /* Receive messages */
399     ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status);
400     if (flag){
401       ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRQ(ierr);
402       proc_id = r_status.MPI_SOURCE;
403       ierr = MPI_Irecv(data2,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr);
404       ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr);
405       if (len > 1+is_max){ /* Add data2 into data */
406         data2_i = data2 + 1 + is_max;
407         for (i=0; i<is_max; i++){
408           table_i = table[i];
409           data_i  = data + 1 + is_max + Mbs*i;
410           isz     = data[1+i];
411           for (j=0; j<data2[1+i]; j++){
412             col = data2_i[j];
413             if (!PetscBTLookupSet(table_i,col)) {data_i[isz++] = col;}
414           }
415           data[1+i] = isz;
416           if (i < is_max - 1) data2_i += data2[1+i];
417         }
418       }
419       k++;
420     }
421   }
422   ierr = PetscFree(data2);CHKERRQ(ierr);
423   ierr = PetscFree2(table,t_p);CHKERRQ(ierr);
424 
425   /* phase 1 sends are complete */
426   ierr = PetscMalloc(size*sizeof(MPI_Status),&s_status);CHKERRQ(ierr);
427   if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);}
428   ierr = PetscFree(data1);CHKERRQ(ierr);
429 
430   /* phase 2 sends are complete */
431   if (nrqr){ierr = MPI_Waitall(nrqr,s_waits2,s_status);CHKERRQ(ierr);}
432   ierr = PetscFree2(s_waits1,s_waits2);CHKERRQ(ierr);
433   ierr = PetscFree(s_status);CHKERRQ(ierr);
434 
435   /* 5. Create new is[] */
436   /*--------------------*/
437   for (i=0; i<is_max; i++) {
438     data_i = data + 1 + is_max + Mbs*i;
439     ierr = ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
440   }
441   for (k=0; k<=nodata2; k++){
442     ierr = PetscFree(odata2_ptr[k]);CHKERRQ(ierr);
443   }
444   ierr = PetscFree(odata2_ptr);CHKERRQ(ierr);
445 
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Local"
451 /*
452    MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do
453        the work on the local processor.
454 
455      Inputs:
456       C      - MAT_MPISBAIJ;
457       data   - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format.
458       whose  - whose is[] to be processed,
459                MINE:  this processor's is[]
460                OTHER: other processor's is[]
461      Output:
462        nidx  - whose = MINE:
463                      holds input and newly found indices in the same format as data
464                whose = OTHER:
465                      only holds the newly found indices
466        table - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'.
467 */
468 /* Would computation be reduced by swapping the loop 'for each is' and 'for each row'? */
469 static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat C,PetscInt *data,PetscInt whose,PetscInt *nidx,PetscBT *table)
470 {
471   Mat_MPISBAIJ   *c = (Mat_MPISBAIJ*)C->data;
472   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)(c->A)->data;
473   Mat_SeqBAIJ    *b = (Mat_SeqBAIJ*)(c->B)->data;
474   PetscErrorCode ierr;
475   PetscInt       row,mbs,Mbs,*nidx_i,col,col_max,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l;
476   PetscInt       a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n;
477   PetscBT        table0;  /* mark the indices of input is[] for look up */
478   PetscBT        table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */
479 
480   PetscFunctionBegin;
481   Mbs = c->Mbs; mbs = a->mbs;
482   ai = a->i; aj = a->j;
483   bi = b->i; bj = b->j;
484   garray = c->garray;
485   rstart = c->rstartbs;
486   is_max = data[0];
487 
488   ierr = PetscBTCreate(Mbs,table0);CHKERRQ(ierr);
489 
490   nidx[0] = is_max;
491   idx_i   = data + is_max + 1; /* ptr to input is[0] array */
492   nidx_i  = nidx + is_max + 1; /* ptr to output is[0] array */
493   for (i=0; i<is_max; i++) { /* for each is */
494     isz  = 0;
495     n = data[1+i]; /* size of input is[i] */
496 
497     /* initialize and set table_i(mark idx and nidx) and table0(only mark idx) */
498     if (whose == MINE){ /* process this processor's is[] */
499       table_i = table[i];
500       nidx_i  = nidx + 1+ is_max + Mbs*i;
501     } else {            /* process other processor's is[] - only use one temp table */
502       table_i = table[0];
503     }
504     ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr);
505     ierr = PetscBTMemzero(Mbs,table0);CHKERRQ(ierr);
506     if (n==0) {
507        nidx[1+i] = 0; /* size of new is[i] */
508        continue;
509     }
510 
511     isz0 = 0; col_max = 0;
512     for (j=0; j<n; j++){
513       col = idx_i[j];
514       if (col >= Mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index col %D >= Mbs %D",col,Mbs);
515       if(!PetscBTLookupSet(table_i,col)) {
516         ierr = PetscBTSet(table0,col);CHKERRQ(ierr);
517         if (whose == MINE) {nidx_i[isz0] = col;}
518         if (col_max < col) col_max = col;
519         isz0++;
520       }
521     }
522 
523     if (whose == MINE) {isz = isz0;}
524     k = 0;  /* no. of indices from input is[i] that have been examined */
525     for (row=0; row<mbs; row++){
526       a_start = ai[row]; a_end = ai[row+1];
527       b_start = bi[row]; b_end = bi[row+1];
528       if (PetscBTLookup(table0,row+rstart)){ /* row is on input is[i]:
529                                                 do row search: collect all col in this row */
530         for (l = a_start; l<a_end ; l++){ /* Amat */
531           col = aj[l] + rstart;
532           if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;}
533         }
534         for (l = b_start; l<b_end ; l++){ /* Bmat */
535           col = garray[bj[l]];
536           if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;}
537         }
538         k++;
539         if (k >= isz0) break; /* for (row=0; row<mbs; row++) */
540       } else { /* row is not on input is[i]:
541                   do col serach: add row onto nidx_i if there is a col in nidx_i */
542         for (l = a_start; l<a_end ; l++){ /* Amat */
543           col = aj[l] + rstart;
544           if (col > col_max) break;
545           if (PetscBTLookup(table0,col)){
546             if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;}
547             break; /* for l = start; l<end ; l++) */
548           }
549         }
550         for (l = b_start; l<b_end ; l++){ /* Bmat */
551           col = garray[bj[l]];
552           if (col > col_max) break;
553           if (PetscBTLookup(table0,col)){
554             if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;}
555             break; /* for l = start; l<end ; l++) */
556           }
557         }
558       }
559     }
560 
561     if (i < is_max - 1){
562       idx_i  += n;   /* ptr to input is[i+1] array */
563       nidx_i += isz; /* ptr to output is[i+1] array */
564     }
565     nidx[1+i] = isz; /* size of new is[i] */
566   } /* for each is */
567   ierr = PetscBTDestroy(table0);CHKERRQ(ierr);
568 
569   PetscFunctionReturn(0);
570 }
571 
572 
573