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