Actual source code: sbaijov.c
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>
9: static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat,PetscInt,IS*);
10: static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat,PetscInt*,PetscInt,PetscInt*,PetscBT*);
12: PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat C,PetscInt is_max,IS is[],PetscInt ov)
13: {
14: PetscInt i,N=C->cmap->N, bs=C->rmap->bs,M=C->rmap->N,Mbs=M/bs,*nidx,isz,iov;
15: IS *is_new,*is_row;
16: Mat *submats;
17: Mat_MPISBAIJ *c=(Mat_MPISBAIJ*)C->data;
18: Mat_SeqSBAIJ *asub_i;
19: PetscBT table;
20: PetscInt *ai,brow,nz,nis,l,nmax,nstages_local,nstages,max_no,pos;
21: const PetscInt *idx;
22: PetscBool flg;
24: PetscMalloc1(is_max,&is_new);
25: /* Convert the indices into block format */
26: ISCompressIndicesGeneral(N,C->rmap->n,bs,is_max,is,is_new);
29: /* ----- previous non-scalable implementation ----- */
30: flg = PETSC_FALSE;
31: PetscOptionsHasName(NULL,NULL, "-IncreaseOverlap_old", &flg);
32: if (flg) { /* previous non-scalable implementation */
33: printf("use previous non-scalable implementation...\n");
34: for (i=0; i<ov; ++i) {
35: MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);
36: }
37: } else { /* implementation using modified BAIJ routines */
39: PetscMalloc1(Mbs+1,&nidx);
40: PetscBTCreate(Mbs,&table); /* for column search */
42: /* Create is_row */
43: PetscMalloc1(is_max,&is_row);
44: ISCreateStride(PETSC_COMM_SELF,Mbs,0,1,&is_row[0]);
46: for (i=1; i<is_max; i++) {
47: is_row[i] = is_row[0]; /* reuse is_row[0] */
48: }
50: /* Allocate memory to hold all the submatrices - Modified from MatCreateSubMatrices_MPIBAIJ() */
51: PetscMalloc1(is_max+1,&submats);
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);
58: /* Make sure every processor loops through the nstages */
59: MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));
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: c->ijonly = PETSC_TRUE; /* only matrix data structures are requested */
68: /* The resulting submatrices should be BAIJ, not SBAIJ, hence we change this value to trigger that */
69: PetscStrcpy(((PetscObject)c->A)->type_name,MATSEQBAIJ);
70: MatCreateSubMatrices_MPIBAIJ_local(C,max_no,is_row+pos,is_new+pos,MAT_INITIAL_MATRIX,submats+pos);
71: PetscStrcpy(((PetscObject)c->A)->type_name,MATSEQSBAIJ);
72: pos += max_no;
73: }
75: /* 2) Row search */
76: MatIncreaseOverlap_MPIBAIJ_Once(C,is_max,is_new);
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;
83: /* put is_new obtained from MatIncreaseOverlap_MPIBAIJ() to table */
84: PetscBTMemzero(Mbs,table);
86: ISGetIndices(is_new[i],&idx);
87: ISGetLocalSize(is_new[i],&nis);
88: for (l=0; l<nis; l++) {
89: PetscBTSet(table,idx[l]);
90: nidx[l] = idx[l];
91: }
92: isz = nis;
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: ISRestoreIndices(is_new[i],&idx);
102: ISDestroy(&is_new[i]);
104: /* create updated is_new */
105: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,is_new+i);
106: }
108: /* Free tmp spaces */
109: for (i=0; i<is_max; i++) {
110: MatDestroy(&submats[i]);
111: }
112: }
114: PetscBTDestroy(&table);
115: PetscFree(submats);
116: ISDestroy(&is_row[0]);
117: PetscFree(is_row);
118: PetscFree(nidx);
120: }
122: for (i=0; i<is_max; i++) ISDestroy(&is[i]);
123: ISExpandIndicesGeneral(N,N,bs,is_max,is_new,is);
125: for (i=0; i<is_max; i++) ISDestroy(&is_new[i]);
126: PetscFree(is_new);
127: return 0;
128: }
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: static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat C,PetscInt is_max,IS is[])
151: {
152: Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data;
153: PetscMPIInt size,rank,tag1,tag2,*len_s,nrqr,nrqs,*id_r1,*len_r1,flag,len,*iwork;
154: const PetscInt *idx_i;
155: PetscInt idx,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i;
156: PetscInt Mbs,i,j,k,*odata1,*odata2;
157: PetscInt proc_id,**odata2_ptr,*ctable=NULL,*btable,len_max,len_est;
158: PetscInt proc_end=0,len_unused,nodata2;
159: PetscInt ois_max; /* max no of is[] in each of processor */
160: char *t_p;
161: MPI_Comm comm;
162: MPI_Request *s_waits1,*s_waits2,r_req;
163: MPI_Status *s_status,r_status;
164: PetscBT *table; /* mark indices of this processor's is[] */
165: PetscBT table_i;
166: PetscBT otable; /* mark indices of other processors' is[] */
167: PetscInt bs=C->rmap->bs,Bn = c->B->cmap->n,Bnbs = Bn/bs,*Bowners;
168: IS garray_local,garray_gl;
170: PetscObjectGetComm((PetscObject)C,&comm);
171: size = c->size;
172: rank = c->rank;
173: Mbs = c->Mbs;
175: PetscObjectGetNewTag((PetscObject)C,&tag1);
176: PetscObjectGetNewTag((PetscObject)C,&tag2);
178: /* create tables used in
179: step 1: table[i] - mark c->garray of proc [i]
180: step 3: table[i] - mark indices of is[i] when whose=MINE
181: table[0] - mark incideces of is[] when whose=OTHER */
182: len = PetscMax(is_max, size);
183: PetscMalloc2(len,&table,(Mbs/PETSC_BITS_PER_BYTE+1)*len,&t_p);
184: for (i=0; i<len; i++) {
185: table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
186: }
188: MPIU_Allreduce(&is_max,&ois_max,1,MPIU_INT,MPI_MAX,comm);
190: /* 1. Send this processor's is[] to other processors */
191: /*---------------------------------------------------*/
192: /* allocate spaces */
193: PetscMalloc1(is_max,&n);
194: len = 0;
195: for (i=0; i<is_max; i++) {
196: ISGetLocalSize(is[i],&n[i]);
197: len += n[i];
198: }
199: if (!len) {
200: is_max = 0;
201: } else {
202: len += 1 + is_max; /* max length of data1 for one processor */
203: }
205: PetscMalloc1(size*len+1,&data1);
206: PetscMalloc1(size,&data1_start);
207: for (i=0; i<size; i++) data1_start[i] = data1 + i*len;
209: PetscMalloc4(size,&len_s,size,&btable,size,&iwork,size+1,&Bowners);
211: /* gather c->garray from all processors */
212: ISCreateGeneral(comm,Bnbs,c->garray,PETSC_COPY_VALUES,&garray_local);
213: ISAllGather(garray_local, &garray_gl);
214: ISDestroy(&garray_local);
215: MPI_Allgather(&Bnbs,1,MPIU_INT,Bowners+1,1,MPIU_INT,comm);
217: Bowners[0] = 0;
218: for (i=0; i<size; i++) Bowners[i+1] += Bowners[i];
220: if (is_max) {
221: /* hash table ctable which maps c->row to proc_id) */
222: PetscMalloc1(Mbs,&ctable);
223: for (proc_id=0,j=0; proc_id<size; proc_id++) {
224: for (; j<C->rmap->range[proc_id+1]/bs; j++) ctable[j] = proc_id;
225: }
227: /* hash tables marking c->garray */
228: ISGetIndices(garray_gl,&idx_i);
229: for (i=0; i<size; i++) {
230: table_i = table[i];
231: PetscBTMemzero(Mbs,table_i);
232: for (j = Bowners[i]; j<Bowners[i+1]; j++) { /* go through B cols of proc[i]*/
233: PetscBTSet(table_i,idx_i[j]);
234: }
235: }
236: ISRestoreIndices(garray_gl,&idx_i);
237: } /* if (is_max) */
238: ISDestroy(&garray_gl);
240: /* evaluate communication - mesg to who, length, and buffer space */
241: for (i=0; i<size; i++) len_s[i] = 0;
243: /* header of data1 */
244: for (proc_id=0; proc_id<size; proc_id++) {
245: iwork[proc_id] = 0;
246: *data1_start[proc_id] = is_max;
247: data1_start[proc_id]++;
248: for (j=0; j<is_max; j++) {
249: if (proc_id == rank) {
250: *data1_start[proc_id] = n[j];
251: } else {
252: *data1_start[proc_id] = 0;
253: }
254: data1_start[proc_id]++;
255: }
256: }
258: for (i=0; i<is_max; i++) {
259: ISGetIndices(is[i],&idx_i);
260: for (j=0; j<n[i]; j++) {
261: idx = idx_i[j];
262: *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */
263: proc_end = ctable[idx];
264: for (proc_id=0; proc_id<=proc_end; proc_id++) { /* for others to process */
265: if (proc_id == rank) continue; /* done before this loop */
266: if (proc_id < proc_end && !PetscBTLookup(table[proc_id],idx)) continue; /* no need for sending idx to [proc_id] */
267: *data1_start[proc_id] = idx; data1_start[proc_id]++;
268: len_s[proc_id]++;
269: }
270: }
271: /* update header data */
272: for (proc_id=0; proc_id<size; proc_id++) {
273: if (proc_id== rank) continue;
274: *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id];
275: iwork[proc_id] = len_s[proc_id];
276: }
277: ISRestoreIndices(is[i],&idx_i);
278: }
280: nrqs = 0; nrqr = 0;
281: for (i=0; i<size; i++) {
282: data1_start[i] = data1 + i*len;
283: if (len_s[i]) {
284: nrqs++;
285: len_s[i] += 1 + is_max; /* add no. of header msg */
286: }
287: }
289: for (i=0; i<is_max; i++) {
290: ISDestroy(&is[i]);
291: }
292: PetscFree(n);
293: PetscFree(ctable);
295: /* Determine the number of messages to expect, their lengths, from from-ids */
296: PetscGatherNumberOfMessages(comm,NULL,len_s,&nrqr);
297: PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);
299: /* Now post the sends */
300: PetscMalloc2(size,&s_waits1,size,&s_waits2);
301: k = 0;
302: for (proc_id=0; proc_id<size; proc_id++) { /* send data1 to processor [proc_id] */
303: if (len_s[proc_id]) {
304: MPI_Isend(data1_start[proc_id],len_s[proc_id],MPIU_INT,proc_id,tag1,comm,s_waits1+k);
305: k++;
306: }
307: }
309: /* 2. Receive other's is[] and process. Then send back */
310: /*-----------------------------------------------------*/
311: len = 0;
312: for (i=0; i<nrqr; i++) {
313: if (len_r1[i] > len) len = len_r1[i];
314: }
315: PetscFree(len_r1);
316: PetscFree(id_r1);
318: for (proc_id=0; proc_id<size; proc_id++) len_s[proc_id] = iwork[proc_id] = 0;
320: PetscMalloc1(len+1,&odata1);
321: PetscMalloc1(size,&odata2_ptr);
322: PetscBTCreate(Mbs,&otable);
324: len_max = ois_max*(Mbs+1); /* max space storing all is[] for each receive */
325: len_est = 2*len_max; /* estimated space of storing is[] for all receiving messages */
326: PetscMalloc1(len_est+1,&odata2);
327: nodata2 = 0; /* nodata2+1: num of PetscMalloc(,&odata2_ptr[]) called */
329: odata2_ptr[nodata2] = odata2;
331: len_unused = len_est; /* unused space in the array odata2_ptr[nodata2]-- needs to be >= len_max */
333: k = 0;
334: while (k < nrqr) {
335: /* Receive messages */
336: MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);
337: if (flag) {
338: MPI_Get_count(&r_status,MPIU_INT,&len);
339: proc_id = r_status.MPI_SOURCE;
340: MPI_Irecv(odata1,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);
341: MPI_Wait(&r_req,&r_status);
343: /* Process messages */
344: /* make sure there is enough unused space in odata2 array */
345: if (len_unused < len_max) { /* allocate more space for odata2 */
346: PetscMalloc1(len_est+1,&odata2);
348: odata2_ptr[++nodata2] = odata2;
350: len_unused = len_est;
351: }
353: MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,odata2,&otable);
354: len = 1 + odata2[0];
355: for (i=0; i<odata2[0]; i++) len += odata2[1 + i];
357: /* Send messages back */
358: MPI_Isend(odata2,len,MPIU_INT,proc_id,tag2,comm,s_waits2+k);
359: k++;
360: odata2 += len;
361: len_unused -= len;
362: len_s[proc_id] = len; /* num of messages sending back to [proc_id] by this proc */
363: }
364: }
365: PetscFree(odata1);
366: PetscBTDestroy(&otable);
368: /* 3. Do local work on this processor's is[] */
369: /*-------------------------------------------*/
370: /* make sure there is enough unused space in odata2(=data) array */
371: len_max = is_max*(Mbs+1); /* max space storing all is[] for this processor */
372: if (len_unused < len_max) { /* allocate more space for odata2 */
373: PetscMalloc1(len_est+1,&odata2);
375: odata2_ptr[++nodata2] = odata2;
376: }
378: data = odata2;
379: MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,data,table);
380: PetscFree(data1_start);
382: /* 4. Receive work done on other processors, then merge */
383: /*------------------------------------------------------*/
384: /* get max number of messages that this processor expects to recv */
385: MPIU_Allreduce(len_s,iwork,size,MPI_INT,MPI_MAX,comm);
386: PetscMalloc1(iwork[rank]+1,&data2);
387: PetscFree4(len_s,btable,iwork,Bowners);
389: k = 0;
390: while (k < nrqs) {
391: /* Receive messages */
392: MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status);
393: if (flag) {
394: MPI_Get_count(&r_status,MPIU_INT,&len);
396: proc_id = r_status.MPI_SOURCE;
398: MPI_Irecv(data2,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);
399: MPI_Wait(&r_req,&r_status);
400: if (len > 1+is_max) { /* Add data2 into data */
401: data2_i = data2 + 1 + is_max;
402: for (i=0; i<is_max; i++) {
403: table_i = table[i];
404: data_i = data + 1 + is_max + Mbs*i;
405: isz = data[1+i];
406: for (j=0; j<data2[1+i]; j++) {
407: col = data2_i[j];
408: if (!PetscBTLookupSet(table_i,col)) data_i[isz++] = col;
409: }
410: data[1+i] = isz;
411: if (i < is_max - 1) data2_i += data2[1+i];
412: }
413: }
414: k++;
415: }
416: }
417: PetscFree(data2);
418: PetscFree2(table,t_p);
420: /* phase 1 sends are complete */
421: PetscMalloc1(size,&s_status);
422: if (nrqs) MPI_Waitall(nrqs,s_waits1,s_status);
423: PetscFree(data1);
425: /* phase 2 sends are complete */
426: if (nrqr) MPI_Waitall(nrqr,s_waits2,s_status);
427: PetscFree2(s_waits1,s_waits2);
428: PetscFree(s_status);
430: /* 5. Create new is[] */
431: /*--------------------*/
432: for (i=0; i<is_max; i++) {
433: data_i = data + 1 + is_max + Mbs*i;
434: ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,PETSC_COPY_VALUES,is+i);
435: }
436: for (k=0; k<=nodata2; k++) {
437: PetscFree(odata2_ptr[k]);
438: }
439: PetscFree(odata2_ptr);
440: return 0;
441: }
443: /*
444: MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do
445: the work on the local processor.
447: Inputs:
448: C - MAT_MPISBAIJ;
449: data - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format.
450: whose - whose is[] to be processed,
451: MINE: this processor's is[]
452: OTHER: other processor's is[]
453: Output:
454: nidx - whose = MINE:
455: holds input and newly found indices in the same format as data
456: whose = OTHER:
457: only holds the newly found indices
458: table - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'.
459: */
460: /* Would computation be reduced by swapping the loop 'for each is' and 'for each row'? */
461: static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat C,PetscInt *data,PetscInt whose,PetscInt *nidx,PetscBT *table)
462: {
463: Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data;
464: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(c->A)->data;
465: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(c->B)->data;
466: PetscInt row,mbs,Mbs,*nidx_i,col,col_max,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l;
467: PetscInt a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n;
468: PetscBT table0; /* mark the indices of input is[] for look up */
469: PetscBT table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */
471: Mbs = c->Mbs; mbs = a->mbs;
472: ai = a->i; aj = a->j;
473: bi = b->i; bj = b->j;
474: garray = c->garray;
475: rstart = c->rstartbs;
476: is_max = data[0];
478: PetscBTCreate(Mbs,&table0);
480: nidx[0] = is_max;
481: idx_i = data + is_max + 1; /* ptr to input is[0] array */
482: nidx_i = nidx + is_max + 1; /* ptr to output is[0] array */
483: for (i=0; i<is_max; i++) { /* for each is */
484: isz = 0;
485: n = data[1+i]; /* size of input is[i] */
487: /* initialize and set table_i(mark idx and nidx) and table0(only mark idx) */
488: if (whose == MINE) { /* process this processor's is[] */
489: table_i = table[i];
490: nidx_i = nidx + 1+ is_max + Mbs*i;
491: } else { /* process other processor's is[] - only use one temp table */
492: table_i = table[0];
493: }
494: PetscBTMemzero(Mbs,table_i);
495: PetscBTMemzero(Mbs,table0);
496: if (n==0) {
497: nidx[1+i] = 0; /* size of new is[i] */
498: continue;
499: }
501: isz0 = 0; col_max = 0;
502: for (j=0; j<n; j++) {
503: col = idx_i[j];
505: if (!PetscBTLookupSet(table_i,col)) {
506: PetscBTSet(table0,col);
507: if (whose == MINE) nidx_i[isz0] = col;
508: if (col_max < col) col_max = col;
509: isz0++;
510: }
511: }
513: if (whose == MINE) isz = isz0;
514: k = 0; /* no. of indices from input is[i] that have been examined */
515: for (row=0; row<mbs; row++) {
516: a_start = ai[row]; a_end = ai[row+1];
517: b_start = bi[row]; b_end = bi[row+1];
518: if (PetscBTLookup(table0,row+rstart)) { /* row is on input is[i]:
519: do row search: collect all col in this row */
520: for (l = a_start; l<a_end ; l++) { /* Amat */
521: col = aj[l] + rstart;
522: if (!PetscBTLookupSet(table_i,col)) nidx_i[isz++] = col;
523: }
524: for (l = b_start; l<b_end ; l++) { /* Bmat */
525: col = garray[bj[l]];
526: if (!PetscBTLookupSet(table_i,col)) nidx_i[isz++] = col;
527: }
528: k++;
529: if (k >= isz0) break; /* for (row=0; row<mbs; row++) */
530: } else { /* row is not on input is[i]:
531: do col serach: add row onto nidx_i if there is a col in nidx_i */
532: for (l = a_start; l<a_end; l++) { /* Amat */
533: col = aj[l] + rstart;
534: if (col > col_max) break;
535: if (PetscBTLookup(table0,col)) {
536: if (!PetscBTLookupSet(table_i,row+rstart)) nidx_i[isz++] = row+rstart;
537: break; /* for l = start; l<end ; l++) */
538: }
539: }
540: for (l = b_start; l<b_end; l++) { /* Bmat */
541: col = garray[bj[l]];
542: if (col > col_max) break;
543: if (PetscBTLookup(table0,col)) {
544: if (!PetscBTLookupSet(table_i,row+rstart)) nidx_i[isz++] = row+rstart;
545: break; /* for l = start; l<end ; l++) */
546: }
547: }
548: }
549: }
551: if (i < is_max - 1) {
552: idx_i += n; /* ptr to input is[i+1] array */
553: nidx_i += isz; /* ptr to output is[i+1] array */
554: }
555: nidx[1+i] = isz; /* size of new is[i] */
556: } /* for each is */
557: PetscBTDestroy(&table0);
558: return 0;
559: }