Actual source code: sbaijov.c

  1: #define PETSCMAT_DLL

  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

 10: static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat,PetscInt,IS*);
 11: static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat,PetscInt*,PetscInt,PetscInt*,PetscBT*);

 15: PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat C,PetscInt is_max,IS is[],PetscInt ov)
 16: {
 18:   PetscInt       i,N=C->N, bs=C->bs;
 19:   IS             *is_new;

 22:   PetscMalloc(is_max*sizeof(IS),&is_new);
 23:   /* Convert the indices into block format */
 24:   ISCompressIndicesGeneral(N,bs,is_max,is,is_new);
 25:   if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
 26:   for (i=0; i<ov; ++i) {
 27:     MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);
 28:   }
 29:   for (i=0; i<is_max; i++) {ISDestroy(is[i]);}
 30:   ISExpandIndicesGeneral(N,bs,is_max,is_new,is);
 31:   for (i=0; i<is_max; i++) {ISDestroy(is_new[i]);}
 32:   PetscFree(is_new);
 33:   return(0);
 34: }

 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: */
 58: static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat C,PetscInt is_max,IS is[])
 59: {
 60:   Mat_MPISBAIJ  *c = (Mat_MPISBAIJ*)C->data;
 62:   PetscMPIInt    size,rank,tag1,tag2,*len_s,nrqr,nrqs,*id_r1,*len_r1,flag,len;
 63:   PetscInt       idx,*idx_i,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i,
 64:                  Mbs,i,j,k,*odata1,*odata2,
 65:                  proc_id,**odata2_ptr,*ctable=0,*btable,len_max,len_est;
 66:   PetscInt       proc_end=0,*iwork,len_unused,nodata2;
 67:   PetscInt       ois_max; /* max no of is[] in each of processor */
 68:   char           *t_p;
 69:   MPI_Comm       comm;
 70:   MPI_Request    *s_waits1,*s_waits2,r_req;
 71:   MPI_Status     *s_status,r_status;
 72:   PetscBT        *table;  /* mark indices of this processor's is[] */
 73:   PetscBT        table_i;
 74:   PetscBT        otable; /* mark indices of other processors' is[] */
 75:   PetscInt       bs=C->bs,Bn = c->B->n,Bnbs = Bn/bs,*Bowners;
 76:   IS             garray_local,garray_gl;

 79:   comm = C->comm;
 80:   size = c->size;
 81:   rank = c->rank;
 82:   Mbs  = c->Mbs;

 84:   PetscObjectGetNewTag((PetscObject)C,&tag1);
 85:   PetscObjectGetNewTag((PetscObject)C,&tag2);

 87:   /* create tables used in
 88:      step 1: table[i] - mark c->garray of proc [i]
 89:      step 3: table[i] - mark indices of is[i] when whose=MINE     
 90:              table[0] - mark incideces of is[] when whose=OTHER */
 91:   len = PetscMax(is_max, size);
 92:   len_max = len*sizeof(PetscBT) + (Mbs/PETSC_BITS_PER_BYTE+1)*len*sizeof(char) + 1;
 93:   PetscMalloc(len_max,&table);
 94:   t_p  = (char *)(table + len);
 95:   for (i=0; i<len; i++) {
 96:     table[i]  = t_p  + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
 97:   }

 99:   MPI_Allreduce(&is_max,&ois_max,1,MPIU_INT,MPI_MAX,comm);
100: 
101:   /* 1. Send this processor's is[] to other processors */
102:   /*---------------------------------------------------*/
103:   /* allocate spaces */
104:   PetscMalloc(is_max*sizeof(PetscInt),&n);
105:   len = 0;
106:   for (i=0; i<is_max; i++) {
107:     ISGetLocalSize(is[i],&n[i]);
108:     len += n[i];
109:   }
110:   if (!len) {
111:     is_max = 0;
112:   } else {
113:     len += 1 + is_max; /* max length of data1 for one processor */
114:   }

116: 
117:   PetscMalloc((size*len+1)*sizeof(PetscInt),&data1);
118:   PetscMalloc(size*sizeof(PetscInt*),&data1_start);
119:   for (i=0; i<size; i++) data1_start[i] = data1 + i*len;

121:   PetscMalloc((size*4+1)*sizeof(PetscInt),&len_s);
122:   btable  = (PetscInt*)(len_s + size);
123:   iwork   = btable + size;
124:   Bowners = iwork + size;

126:   /* gather c->garray from all processors */
127:   ISCreateGeneral(comm,Bnbs,c->garray,&garray_local);
128:   ISAllGather(garray_local, &garray_gl);
129:   ISDestroy(garray_local);
130:   MPI_Allgather(&Bnbs,1,MPIU_INT,Bowners+1,1,MPIU_INT,comm);
131:   Bowners[0] = 0;
132:   for (i=0; i<size; i++) Bowners[i+1] += Bowners[i];
133: 
134:   if (is_max){
135:     /* hash table ctable which maps c->row to proc_id) */
136:     PetscMalloc(Mbs*sizeof(PetscInt),&ctable);
137:     for (proc_id=0,j=0; proc_id<size; proc_id++) {
138:       for (; j<c->rowners[proc_id+1]; j++) {
139:         ctable[j] = proc_id;
140:       }
141:     }

143:     /* hash tables marking c->garray */
144:     ISGetIndices(garray_gl,&idx_i);
145:     for (i=0; i<size; i++){
146:       table_i = table[i];
147:       PetscBTMemzero(Mbs,table_i);
148:       for (j = Bowners[i]; j<Bowners[i+1]; j++){ /* go through B cols of proc[i]*/
149:         PetscBTSet(table_i,idx_i[j]);
150:       }
151:     }
152:     ISRestoreIndices(garray_gl,&idx_i);
153:   }  /* if (is_max) */
154:   ISDestroy(garray_gl);

156:   /* evaluate communication - mesg to who, length, and buffer space */
157:   for (i=0; i<size; i++) len_s[i] = 0;
158: 
159:   /* header of data1 */
160:   for (proc_id=0; proc_id<size; proc_id++){
161:     iwork[proc_id] = 0;
162:     *data1_start[proc_id] = is_max;
163:     data1_start[proc_id]++;
164:     for (j=0; j<is_max; j++) {
165:       if (proc_id == rank){
166:         *data1_start[proc_id] = n[j];
167:       } else {
168:         *data1_start[proc_id] = 0;
169:       }
170:       data1_start[proc_id]++;
171:     }
172:   }
173: 
174:   for (i=0; i<is_max; i++) {
175:     ISGetIndices(is[i],&idx_i);
176:     for (j=0; j<n[i]; j++){
177:       idx = idx_i[j];
178:       *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */
179:       proc_end = ctable[idx];
180:       for (proc_id=0;  proc_id<=proc_end; proc_id++){ /* for others to process */
181:         if (proc_id == rank ) continue; /* done before this loop */
182:         if (proc_id < proc_end && !PetscBTLookup(table[proc_id],idx))
183:           continue;   /* no need for sending idx to [proc_id] */
184:         *data1_start[proc_id] = idx; data1_start[proc_id]++;
185:         len_s[proc_id]++;
186:       }
187:     }
188:     /* update header data */
189:     for (proc_id=0; proc_id<size; proc_id++){
190:       if (proc_id== rank) continue;
191:       *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id];
192:       iwork[proc_id] = len_s[proc_id] ;
193:     }
194:     ISRestoreIndices(is[i],&idx_i);
195:   }

197:   nrqs = 0; nrqr = 0;
198:   for (i=0; i<size; i++){
199:     data1_start[i] = data1 + i*len;
200:     if (len_s[i]){
201:       nrqs++;
202:       len_s[i] += 1 + is_max; /* add no. of header msg */
203:     }
204:   }

206:   for (i=0; i<is_max; i++) {
207:     ISDestroy(is[i]);
208:   }
209:   PetscFree(n);
210:   if (ctable){PetscFree(ctable);}

212:   /* Determine the number of messages to expect, their lengths, from from-ids */
213:   PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&nrqr);
214:   PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);
215: 
216:   /*  Now  post the sends */
217:   PetscMalloc(2*size*sizeof(MPI_Request),&s_waits1);
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:       MPI_Isend(data1_start[proc_id],len_s[proc_id],MPIU_INT,proc_id,tag1,comm,s_waits1+k);
223:       k++;
224:     }
225:   }

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:   }
233:   PetscFree(len_r1);
234:   PetscFree(id_r1);

236:   for (proc_id=0; proc_id<size; proc_id++)
237:     len_s[proc_id] = iwork[proc_id] = 0;
238: 
239:   PetscMalloc((len+1)*sizeof(PetscInt),&odata1);
240:   PetscMalloc(size*sizeof(PetscInt**),&odata2_ptr);
241:   PetscBTCreate(Mbs,otable);

243:   len_max = ois_max*(Mbs+1);  /* max space storing all is[] for each receive */
244:   len_est = 2*len_max; /* estimated space of storing is[] for all receiving messages */
245:   PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);
246:   nodata2 = 0;       /* nodata2+1: num of PetscMalloc(,&odata2_ptr[]) called */
247:   odata2_ptr[nodata2] = odata2;
248:   len_unused = len_est; /* unused space in the array odata2_ptr[nodata2]-- needs to be >= len_max  */
249: 
250:   k = 0;
251:   while (k < nrqr){
252:     /* Receive messages */
253:     MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);
254:     if (flag){
255:       MPI_Get_count(&r_status,MPIU_INT,&len);
256:       proc_id = r_status.MPI_SOURCE;
257:       MPI_Irecv(odata1,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);
258:       MPI_Wait(&r_req,&r_status);

260:       /*  Process messages */
261:       /*  make sure there is enough unused space in odata2 array */
262:       if (len_unused < len_max){ /* allocate more space for odata2 */
263:         PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);
264:         odata2_ptr[++nodata2] = odata2;
265:         len_unused = len_est;
266:       }

268:       MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,odata2,&otable);
269:       len = 1 + odata2[0];
270:       for (i=0; i<odata2[0]; i++){
271:         len += odata2[1 + i];
272:       }

274:       /* Send messages back */
275:       MPI_Isend(odata2,len,MPIU_INT,proc_id,tag2,comm,s_waits2+k);
276:       k++;
277:       odata2     += len;
278:       len_unused -= len;
279:       len_s[proc_id] = len; /* num of messages sending back to [proc_id] by this proc */
280:     }
281:   }
282:   PetscFree(odata1);
283:   PetscBTDestroy(otable);

285:   /* 3. Do local work on this processor's is[] */
286:   /*-------------------------------------------*/
287:   /* make sure there is enough unused space in odata2(=data) array */
288:   len_max = is_max*(Mbs+1); /* max space storing all is[] for this processor */
289:   if (len_unused < len_max){ /* allocate more space for odata2 */
290:     PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);
291:     odata2_ptr[++nodata2] = odata2;
292:     len_unused = len_est;
293:   }

295:   data = odata2;
296:   MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,data,table);
297:   PetscFree(data1_start);

299:   /* 4. Receive work done on other processors, then merge */
300:   /*------------------------------------------------------*/
301:   /* get max number of messages that this processor expects to recv */
302:   MPI_Allreduce(len_s,iwork,size,MPIU_INT,MPI_MAX,comm);
303:   PetscMalloc((iwork[rank]+1)*sizeof(PetscInt),&data2);
304:   PetscFree(len_s);

306:   k = 0;
307:   while (k < nrqs){
308:     /* Receive messages */
309:     MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status);
310:     if (flag){
311:       MPI_Get_count(&r_status,MPIU_INT,&len);
312:       proc_id = r_status.MPI_SOURCE;
313:       MPI_Irecv(data2,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);
314:       MPI_Wait(&r_req,&r_status);
315:       if (len > 1+is_max){ /* Add data2 into data */
316:         data2_i = data2 + 1 + is_max;
317:         for (i=0; i<is_max; i++){
318:           table_i = table[i];
319:           data_i  = data + 1 + is_max + Mbs*i;
320:           isz     = data[1+i];
321:           for (j=0; j<data2[1+i]; j++){
322:             col = data2_i[j];
323:             if (!PetscBTLookupSet(table_i,col)) {data_i[isz++] = col;}
324:           }
325:           data[1+i] = isz;
326:           if (i < is_max - 1) data2_i += data2[1+i];
327:         }
328:       }
329:       k++;
330:     }
331:   }
332:   PetscFree(data2);
333:   PetscFree(table);

335:   /* phase 1 sends are complete */
336:   PetscMalloc(size*sizeof(MPI_Status),&s_status);
337:   if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status);}
338:   PetscFree(data1);
339: 
340:   /* phase 2 sends are complete */
341:   if (nrqr){MPI_Waitall(nrqr,s_waits2,s_status);}
342:   PetscFree(s_waits1);
343:   PetscFree(s_status);

345:   /* 5. Create new is[] */
346:   /*--------------------*/
347:   for (i=0; i<is_max; i++) {
348:     data_i = data + 1 + is_max + Mbs*i;
349:     ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,is+i);
350:   }
351:   for (k=0; k<=nodata2; k++){
352:     PetscFree(odata2_ptr[k]);
353:   }
354:   PetscFree(odata2_ptr);

356:   return(0);
357: }

361: /*  
362:    MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do 
363:        the work on the local processor.

365:      Inputs:
366:       C      - MAT_MPISBAIJ;
367:       data   - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format.
368:       whose  - whose is[] to be processed, 
369:                MINE:  this processor's is[]
370:                OTHER: other processor's is[]
371:      Output:  
372:        nidx  - whose = MINE:
373:                      holds input and newly found indices in the same format as data
374:                whose = OTHER:
375:                      only holds the newly found indices
376:        table - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'.
377: */
378: /* Would computation be reduced by swapping the loop 'for each is' and 'for each row'? */
379: static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat C,PetscInt *data,PetscInt whose,PetscInt *nidx,PetscBT *table)
380: {
381:   Mat_MPISBAIJ   *c = (Mat_MPISBAIJ*)C->data;
382:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)(c->A)->data;
383:   Mat_SeqBAIJ    *b = (Mat_SeqBAIJ*)(c->B)->data;
385:   PetscInt       row,mbs,Mbs,*nidx_i,col,col_max,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l;
386:   PetscInt       a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n;
387:   PetscBT        table0;  /* mark the indices of input is[] for look up */
388:   PetscBT        table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */
389: 
391:   Mbs = c->Mbs; mbs = a->mbs;
392:   ai = a->i; aj = a->j;
393:   bi = b->i; bj = b->j;
394:   garray = c->garray;
395:   rstart = c->rstart;
396:   is_max = data[0];

398:   PetscBTCreate(Mbs,table0);
399: 
400:   nidx[0] = is_max;
401:   idx_i   = data + is_max + 1; /* ptr to input is[0] array */
402:   nidx_i  = nidx + is_max + 1; /* ptr to output is[0] array */
403:   for (i=0; i<is_max; i++) { /* for each is */
404:     isz  = 0;
405:     n = data[1+i]; /* size of input is[i] */

407:     /* initialize and set table_i(mark idx and nidx) and table0(only mark idx) */
408:     if (whose == MINE){ /* process this processor's is[] */
409:       table_i = table[i];
410:       nidx_i  = nidx + 1+ is_max + Mbs*i;
411:     } else {            /* process other processor's is[] - only use one temp table */
412:       table_i = table[0];
413:     }
414:     PetscBTMemzero(Mbs,table_i);
415:     PetscBTMemzero(Mbs,table0);
416:     if (n==0) {
417:        nidx[1+i] = 0; /* size of new is[i] */
418:        continue;
419:     }

421:     isz0 = 0; col_max = 0;
422:     for (j=0; j<n; j++){
423:       col = idx_i[j];
424:       if (col >= Mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"index col %D >= Mbs %D",col,Mbs);
425:       if(!PetscBTLookupSet(table_i,col)) {
426:         PetscBTSet(table0,col);
427:         if (whose == MINE) {nidx_i[isz0] = col;}
428:         if (col_max < col) col_max = col;
429:         isz0++;
430:       }
431:     }
432: 
433:     if (whose == MINE) {isz = isz0;}
434:     k = 0;  /* no. of indices from input is[i] that have been examined */
435:     for (row=0; row<mbs; row++){
436:       a_start = ai[row]; a_end = ai[row+1];
437:       b_start = bi[row]; b_end = bi[row+1];
438:       if (PetscBTLookup(table0,row+rstart)){ /* row is on input is[i]:
439:                                                 do row search: collect all col in this row */
440:         for (l = a_start; l<a_end ; l++){ /* Amat */
441:           col = aj[l] + rstart;
442:           if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;}
443:         }
444:         for (l = b_start; l<b_end ; l++){ /* Bmat */
445:           col = garray[bj[l]];
446:           if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;}
447:         }
448:         k++;
449:         if (k >= isz0) break; /* for (row=0; row<mbs; row++) */
450:       } else { /* row is not on input is[i]:
451:                   do col serach: add row onto nidx_i if there is a col in nidx_i */
452:         for (l = a_start; l<a_end ; l++){ /* Amat */
453:           col = aj[l] + rstart;
454:           if (col > col_max) break;
455:           if (PetscBTLookup(table0,col)){
456:             if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;}
457:             break; /* for l = start; l<end ; l++) */
458:           }
459:         }
460:         for (l = b_start; l<b_end ; l++){ /* Bmat */
461:           col = garray[bj[l]];
462:           if (col > col_max) break;
463:           if (PetscBTLookup(table0,col)){
464:             if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;}
465:             break; /* for l = start; l<end ; l++) */
466:           }
467:         }
468:       }
469:     }
470: 
471:     if (i < is_max - 1){
472:       idx_i  += n;   /* ptr to input is[i+1] array */
473:       nidx_i += isz; /* ptr to output is[i+1] array */
474:     }
475:     nidx[1+i] = isz; /* size of new is[i] */
476:   } /* for each is */
477:   PetscBTDestroy(table0);
478: 
479:   return(0);
480: }