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