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