xref: /petsc/src/mat/impls/sbaij/mpi/sbaijov.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1632d0f97SHong Zhang 
2632d0f97SHong Zhang /*
3632d0f97SHong Zhang    Routines to compute overlapping regions of a parallel MPI matrix.
4632d0f97SHong Zhang    Used for finding submatrices that were shared across processors.
5632d0f97SHong Zhang */
6632d0f97SHong Zhang #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
7632d0f97SHong Zhang #include "petscbt.h"
8632d0f97SHong Zhang 
9*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat,int,IS*);
10*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat,int*,int,int*,PetscBT*);
11632d0f97SHong Zhang 
12632d0f97SHong Zhang #undef __FUNCT__
13632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ"
14dfbe8321SBarry Smith PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat C,int is_max,IS is[],int ov)
15632d0f97SHong Zhang {
16d9489beaSHong Zhang   Mat_MPISBAIJ  *c = (Mat_MPISBAIJ*)C->data;
17*6849ba73SBarry Smith   PetscErrorCode ierr;
18*6849ba73SBarry Smith   int           i,N=C->N, bs=c->bs;
19632d0f97SHong Zhang   IS            *is_new;
20632d0f97SHong Zhang 
21632d0f97SHong Zhang   PetscFunctionBegin;
22c910923dSHong Zhang   ierr = PetscMalloc(is_max*sizeof(IS),&is_new);CHKERRQ(ierr);
23632d0f97SHong Zhang   /* Convert the indices into block format */
2427f478b1SHong Zhang   ierr = ISCompressIndicesGeneral(N,bs,is_max,is,is_new);CHKERRQ(ierr);
25632d0f97SHong Zhang   if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");}
26632d0f97SHong Zhang   for (i=0; i<ov; ++i) {
27c910923dSHong Zhang     ierr = MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr);
28632d0f97SHong Zhang   }
29c910923dSHong Zhang   for (i=0; i<is_max; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);}
3027f478b1SHong Zhang   ierr = ISExpandIndicesGeneral(N,bs,is_max,is_new,is);CHKERRQ(ierr);
31c910923dSHong Zhang   for (i=0; i<is_max; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);}
32632d0f97SHong Zhang   ierr = PetscFree(is_new);CHKERRQ(ierr);
33632d0f97SHong Zhang   PetscFunctionReturn(0);
34632d0f97SHong Zhang }
35632d0f97SHong Zhang 
364a69c536SBarry Smith typedef enum {MINE,OTHER} WhoseOwner;
370472cc68SHong Zhang /*  data1, odata1 and odata2 are packed in the format (for communication):
38a2a9f22aSHong Zhang        data[0]          = is_max, no of is
39a2a9f22aSHong Zhang        data[1]          = size of is[0]
40a2a9f22aSHong Zhang         ...
41a2a9f22aSHong Zhang        data[is_max]     = size of is[is_max-1]
42a2a9f22aSHong Zhang        data[is_max + 1] = data(is[0])
43a2a9f22aSHong Zhang         ...
44a2a9f22aSHong Zhang        data[is_max+1+sum(size of is[k]), k=0,...,i-1] = data(is[i])
45a2a9f22aSHong Zhang         ...
460472cc68SHong Zhang    data2 is packed in the format (for creating output is[]):
470472cc68SHong Zhang        data[0]          = is_max, no of is
480472cc68SHong Zhang        data[1]          = size of is[0]
490472cc68SHong Zhang         ...
500472cc68SHong Zhang        data[is_max]     = size of is[is_max-1]
510472cc68SHong Zhang        data[is_max + 1] = data(is[0])
520472cc68SHong Zhang         ...
530472cc68SHong Zhang        data[is_max + 1 + Mbs*i) = data(is[i])
540472cc68SHong Zhang         ...
55a2a9f22aSHong Zhang */
56632d0f97SHong Zhang #undef __FUNCT__
57632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Once"
58*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat C,int is_max,IS is[])
59632d0f97SHong Zhang {
60632d0f97SHong Zhang   Mat_MPISBAIJ  *c = (Mat_MPISBAIJ*)C->data;
61*6849ba73SBarry Smith   PetscErrorCode ierr;
625483b11dSHong Zhang   int         len,idx,*idx_i,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i,
63*6849ba73SBarry Smith               size,rank,Mbs,i,j,k,nrqs,nrqr,*odata1,*odata2,
6410eea884SHong Zhang               tag1,tag2,flag,proc_id,**odata2_ptr,*ctable=0,*btable,len_max,len_est;
6510eea884SHong Zhang   int         *id_r1,*len_r1,proc_end=0,*iwork,*len_s,len_unused,nodata2;
6610eea884SHong Zhang   int         ois_max; /* max no of is[] in each of processor */
67bfc6803cSHong Zhang   char        *t_p;
68632d0f97SHong Zhang   MPI_Comm    comm;
69e8527bf2SHong Zhang   MPI_Request *s_waits1,*s_waits2,r_req;
70632d0f97SHong Zhang   MPI_Status  *s_status,r_status;
7145f43ab7SHong Zhang   PetscBT     *table;  /* mark indices of this processor's is[] */
72fc70d14eSHong Zhang   PetscBT     table_i;
73bfc6803cSHong Zhang   PetscBT     otable; /* mark indices of other processors' is[] */
74e8527bf2SHong Zhang   int         bs=c->bs,Bn = c->B->n,Bnbs = Bn/bs,*Bowners;
7510eea884SHong Zhang   IS          garray_local,garray_gl;
765483b11dSHong Zhang 
77632d0f97SHong Zhang   PetscFunctionBegin;
78632d0f97SHong Zhang 
79632d0f97SHong Zhang   comm = C->comm;
80632d0f97SHong Zhang   size = c->size;
81632d0f97SHong Zhang   rank = c->rank;
82632d0f97SHong Zhang   Mbs  = c->Mbs;
83632d0f97SHong Zhang 
84c910923dSHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
85c910923dSHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr);
86c910923dSHong Zhang 
87430a0127SHong Zhang   /* create tables used in
88430a0127SHong Zhang      step 1: table[i] - mark c->garray of proc [i]
8945f43ab7SHong Zhang      step 3: table[i] - mark indices of is[i] when whose=MINE
90430a0127SHong Zhang              table[0] - mark incideces of is[] when whose=OTHER */
91430a0127SHong Zhang   len = PetscMax(is_max, size);CHKERRQ(ierr);
92430a0127SHong Zhang   len_max = len*sizeof(PetscBT) + (Mbs/PETSC_BITS_PER_BYTE+1)*len*sizeof(char) + 1;
93430a0127SHong Zhang   ierr = PetscMalloc(len_max,&table);CHKERRQ(ierr);
94430a0127SHong Zhang   t_p  = (char *)(table + len);
95430a0127SHong Zhang   for (i=0; i<len; i++) {
96430a0127SHong Zhang     table[i]  = t_p  + (Mbs/PETSC_BITS_PER_BYTE+1)*i;
97430a0127SHong Zhang   }
98430a0127SHong Zhang 
9910eea884SHong Zhang   ierr = MPI_Allreduce(&is_max,&ois_max,1,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
10010eea884SHong Zhang 
1015483b11dSHong Zhang   /* 1. Send this processor's is[] to other processors */
1025483b11dSHong Zhang   /*---------------------------------------------------*/
103e8527bf2SHong Zhang   /* allocate spaces */
1045483b11dSHong Zhang   ierr = PetscMalloc(is_max*sizeof(int),&n);CHKERRQ(ierr);
1055483b11dSHong Zhang   len = 0;
106c910923dSHong Zhang   for (i=0; i<is_max; i++) {
107632d0f97SHong Zhang     ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr);
108632d0f97SHong Zhang     len += n[i];
109632d0f97SHong Zhang   }
110958c9bccSBarry Smith   if (!len) {
1115483b11dSHong Zhang     is_max = 0;
1125483b11dSHong Zhang   } else {
1135483b11dSHong Zhang     len += 1 + is_max; /* max length of data1 for one processor */
1145483b11dSHong Zhang   }
1155483b11dSHong Zhang 
116430a0127SHong Zhang 
1175483b11dSHong Zhang   ierr = PetscMalloc((size*len+1)*sizeof(int),&data1);CHKERRQ(ierr);
1185483b11dSHong Zhang   ierr = PetscMalloc(size*sizeof(int*),&data1_start);CHKERRQ(ierr);
1195483b11dSHong Zhang   for (i=0; i<size; i++) data1_start[i] = data1 + i*len;
1205483b11dSHong Zhang 
12110eea884SHong Zhang   ierr = PetscMalloc((size*4+1)*sizeof(int),&len_s);CHKERRQ(ierr);
122e8527bf2SHong Zhang   btable  = len_s + size;
123e8527bf2SHong Zhang   iwork   = btable + size;
12410eea884SHong Zhang   Bowners = iwork + size;
125e8527bf2SHong Zhang 
12676f244e2SHong Zhang   /* gather c->garray from all processors */
12776f244e2SHong Zhang   ierr = ISCreateGeneral(comm,Bnbs,c->garray,&garray_local);CHKERRQ(ierr);
12876f244e2SHong Zhang   ierr = ISAllGather(garray_local, &garray_gl);CHKERRQ(ierr);
12976f244e2SHong Zhang   ierr = ISDestroy(garray_local);CHKERRQ(ierr);
13076f244e2SHong Zhang   ierr = MPI_Allgather(&Bnbs,1,MPI_INT,Bowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
13176f244e2SHong Zhang   Bowners[0] = 0;
13276f244e2SHong Zhang   for (i=0; i<size; i++) Bowners[i+1] += Bowners[i];
13376f244e2SHong Zhang 
1345483b11dSHong Zhang   if (is_max){
135430a0127SHong Zhang     /* hash table ctable which maps c->row to proc_id) */
1365483b11dSHong Zhang     ierr = PetscMalloc(Mbs*sizeof(int),&ctable);CHKERRQ(ierr);
1375483b11dSHong Zhang     for (proc_id=0,j=0; proc_id<size; proc_id++) {
1385483b11dSHong Zhang       for (; j<c->rowners[proc_id+1]; j++) {
1395483b11dSHong Zhang         ctable[j] = proc_id;
1405483b11dSHong Zhang       }
1415483b11dSHong Zhang     }
1425483b11dSHong Zhang 
143430a0127SHong Zhang     /* hash tables marking c->garray */
14410eea884SHong Zhang     ierr = ISGetIndices(garray_gl,&idx_i);
1455483b11dSHong Zhang     for (i=0; i<size; i++){
146430a0127SHong Zhang       table_i = table[i];
147430a0127SHong Zhang       ierr    = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr);
148430a0127SHong Zhang       for (j = Bowners[i]; j<Bowners[i+1]; j++){ /* go through B cols of proc[i]*/
149430a0127SHong Zhang         ierr = PetscBTSet(table_i,idx_i[j]);CHKERRQ(ierr);
1505483b11dSHong Zhang       }
1515483b11dSHong Zhang     }
15210eea884SHong Zhang     ierr = ISRestoreIndices(garray_gl,&idx_i);CHKERRQ(ierr);
1535483b11dSHong Zhang   }  /* if (is_max) */
15476f244e2SHong Zhang   ierr = ISDestroy(garray_gl);CHKERRQ(ierr);
1555483b11dSHong Zhang 
1565483b11dSHong Zhang   /* evaluate communication - mesg to who, length, and buffer space */
157e8527bf2SHong Zhang   for (i=0; i<size; i++) len_s[i] = 0;
1585483b11dSHong Zhang 
1595483b11dSHong Zhang   /* header of data1 */
1605483b11dSHong Zhang   for (proc_id=0; proc_id<size; proc_id++){
1615483b11dSHong Zhang     iwork[proc_id] = 0;
1625483b11dSHong Zhang     *data1_start[proc_id] = is_max;
1635483b11dSHong Zhang     data1_start[proc_id]++;
1645483b11dSHong Zhang     for (j=0; j<is_max; j++) {
1655483b11dSHong Zhang       if (proc_id == rank){
1665483b11dSHong Zhang         *data1_start[proc_id] = n[j];
1675483b11dSHong Zhang       } else {
1685483b11dSHong Zhang         *data1_start[proc_id] = 0;
1695483b11dSHong Zhang       }
1705483b11dSHong Zhang       data1_start[proc_id]++;
1715483b11dSHong Zhang     }
1725483b11dSHong Zhang   }
1735483b11dSHong Zhang 
1745483b11dSHong Zhang   for (i=0; i<is_max; i++) {
1755483b11dSHong Zhang     ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr);
1765483b11dSHong Zhang     for (j=0; j<n[i]; j++){
1775483b11dSHong Zhang       idx = idx_i[j];
1785483b11dSHong Zhang       *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */
1795483b11dSHong Zhang       proc_end = ctable[idx];
1805483b11dSHong Zhang       for (proc_id=0;  proc_id<=proc_end; proc_id++){ /* for others to process */
181e8527bf2SHong Zhang         if (proc_id == rank ) continue; /* done before this loop */
182430a0127SHong Zhang         if (proc_id < proc_end && !PetscBTLookup(table[proc_id],idx))
183430a0127SHong Zhang           continue;   /* no need for sending idx to [proc_id] */
1845483b11dSHong Zhang         *data1_start[proc_id] = idx; data1_start[proc_id]++;
1855483b11dSHong Zhang         len_s[proc_id]++;
1865483b11dSHong Zhang       }
1875483b11dSHong Zhang     }
1885483b11dSHong Zhang     /* update header data */
1892cfbe0a4SHong Zhang     for (proc_id=0; proc_id<size; proc_id++){
1905483b11dSHong Zhang       if (proc_id== rank) continue;
1915483b11dSHong Zhang       *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id];
1925483b11dSHong Zhang       iwork[proc_id] = len_s[proc_id] ;
1935483b11dSHong Zhang     }
1945483b11dSHong Zhang     ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr);
195e8527bf2SHong Zhang   }
1965483b11dSHong Zhang 
1975483b11dSHong Zhang   nrqs = 0; nrqr = 0;
1985483b11dSHong Zhang   for (i=0; i<size; i++){
1995483b11dSHong Zhang     data1_start[i] = data1 + i*len;
2005483b11dSHong Zhang     if (len_s[i]){
2015483b11dSHong Zhang       nrqs++;
2025483b11dSHong Zhang       len_s[i] += 1 + is_max; /* add no. of header msg */
2035483b11dSHong Zhang     }
2045483b11dSHong Zhang   }
2055483b11dSHong Zhang 
2065483b11dSHong Zhang   for (i=0; i<is_max; i++) {
207a2a9f22aSHong Zhang     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
208632d0f97SHong Zhang   }
209bfc6803cSHong Zhang   ierr = PetscFree(n);CHKERRQ(ierr);
2105483b11dSHong Zhang   if (ctable){ierr = PetscFree(ctable);CHKERRQ(ierr);}
2115483b11dSHong Zhang 
2125483b11dSHong Zhang   /* Determine the number of messages to expect, their lengths, from from-ids */
2135483b11dSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&nrqr);CHKERRQ(ierr);
2145483b11dSHong Zhang   ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);CHKERRQ(ierr);
2155483b11dSHong Zhang   /* ierr = PetscPrintf(PETSC_COMM_SELF, "[%d] nrqs: %d, nrqr: %d\n",rank,nrqs,nrqr); */
216632d0f97SHong Zhang 
217632d0f97SHong Zhang   /*  Now  post the sends */
2185483b11dSHong Zhang   ierr = PetscMalloc(2*size*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
2195483b11dSHong Zhang   s_waits2 = s_waits1 + size;
220632d0f97SHong Zhang   k = 0;
2215483b11dSHong Zhang   for (proc_id=0; proc_id<size; proc_id++){  /* send data1 to processor [proc_id] */
2225483b11dSHong Zhang     if (len_s[proc_id]){
2235483b11dSHong Zhang       ierr = MPI_Isend(data1_start[proc_id],len_s[proc_id],MPI_INT,proc_id,tag1,comm,s_waits1+k);CHKERRQ(ierr);
224632d0f97SHong Zhang       k++;
225632d0f97SHong Zhang     }
226632d0f97SHong Zhang   }
227632d0f97SHong Zhang 
22845f43ab7SHong Zhang   /* 2. Receive other's is[] and process. Then send back */
229bfc6803cSHong Zhang   /*-----------------------------------------------------*/
2305483b11dSHong Zhang   len = 0;
2315483b11dSHong Zhang   for (i=0; i<nrqr; i++){
2325483b11dSHong Zhang     if (len_r1[i] > len)len = len_r1[i];
2335483b11dSHong Zhang     /* ierr = PetscPrintf(PETSC_COMM_SELF, "[%d] expect to recv len=%d from [%d]\n",rank,len_r1[i],id_r1[i]); */
2345483b11dSHong Zhang   }
23545f43ab7SHong Zhang   ierr = PetscFree(len_r1);CHKERRQ(ierr);
23645f43ab7SHong Zhang   ierr = PetscFree(id_r1);CHKERRQ(ierr);
23745f43ab7SHong Zhang 
23845f43ab7SHong Zhang   for (proc_id=0; proc_id<size; proc_id++)
23945f43ab7SHong Zhang     len_s[proc_id] = iwork[proc_id] = 0;
24045f43ab7SHong Zhang 
2415483b11dSHong Zhang   ierr = PetscMalloc((len+1)*sizeof(int),&odata1);CHKERRQ(ierr);
242240e5896SHong Zhang   ierr = PetscMalloc(size*sizeof(int**),&odata2_ptr);CHKERRQ(ierr);
24376f244e2SHong Zhang   ierr = PetscBTCreate(Mbs,otable);CHKERRQ(ierr);
24410eea884SHong Zhang 
24510eea884SHong Zhang   len_max = ois_max*(Mbs+1);  /* max space storing all is[] for each receive */
246240e5896SHong Zhang   len_est = 2*len_max; /* estimated space of storing is[] for all receiving messages */
247240e5896SHong Zhang   ierr = PetscMalloc((len_est+1)*sizeof(int),&odata2);CHKERRQ(ierr);
24810eea884SHong Zhang   nodata2 = 0;       /* nodata2+1: num of PetscMalloc(,&odata2_ptr[]) called */
249240e5896SHong Zhang   odata2_ptr[nodata2] = odata2;
25010eea884SHong Zhang   len_unused = len_est; /* unused space in the array odata2_ptr[nodata2]-- needs to be >= len_max  */
25110eea884SHong Zhang 
252632d0f97SHong Zhang   k = 0;
2535483b11dSHong Zhang   while (k < nrqr){
254632d0f97SHong Zhang     /* Receive messages */
255bfc6803cSHong Zhang     ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);CHKERRQ(ierr);
256632d0f97SHong Zhang     if (flag){
257bfc6803cSHong Zhang       ierr = MPI_Get_count(&r_status,MPI_INT,&len);CHKERRQ(ierr);
258632d0f97SHong Zhang       proc_id = r_status.MPI_SOURCE;
259bfc6803cSHong Zhang       ierr = MPI_Irecv(odata1,len,MPI_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr);
2605483b11dSHong Zhang       ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr);
2615483b11dSHong Zhang       /* ierr = PetscPrintf(PETSC_COMM_SELF, " [%d] recv %d from [%d]\n",rank,len,proc_id); */
262632d0f97SHong Zhang 
263fc70d14eSHong Zhang       /*  Process messages */
264240e5896SHong Zhang       /*  make sure there is enough unused space in odata2 array */
26510eea884SHong Zhang       if (len_unused < len_max){ /* allocate more space for odata2 */
266240e5896SHong Zhang         ierr = PetscMalloc((len_est+1)*sizeof(int),&odata2);CHKERRQ(ierr);
267240e5896SHong Zhang         odata2_ptr[++nodata2] = odata2;
26810eea884SHong Zhang         len_unused = len_est;
269240e5896SHong Zhang         /* ierr = PetscPrintf(PETSC_COMM_SELF, " [%d] 2. Malloc odata2, nodata2: %d\n",rank,nodata2); */
27010eea884SHong Zhang       }
27110eea884SHong Zhang 
27210eea884SHong Zhang       ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,odata2,&otable);CHKERRQ(ierr);
273a2a9f22aSHong Zhang       len = 1 + odata2[0];
274a2a9f22aSHong Zhang       for (i=0; i<odata2[0]; i++){
275a2a9f22aSHong Zhang         len += odata2[1 + i];
276fc70d14eSHong Zhang       }
277632d0f97SHong Zhang 
278632d0f97SHong Zhang       /* Send messages back */
2795483b11dSHong Zhang       ierr = MPI_Isend(odata2,len,MPI_INT,proc_id,tag2,comm,s_waits2+k);CHKERRQ(ierr);
2805483b11dSHong Zhang       /* ierr = PetscPrintf(PETSC_COMM_SELF," [%d] send %d back to [%d] \n",rank,len,proc_id); */
281fc70d14eSHong Zhang       k++;
28210eea884SHong Zhang       odata2     += len;
28310eea884SHong Zhang       len_unused -= len;
28445f43ab7SHong Zhang       len_s[proc_id] = len; /* num of messages sending back to [proc_id] by this proc */
285632d0f97SHong Zhang     }
2865483b11dSHong Zhang   }
2875483b11dSHong Zhang   ierr = PetscFree(odata1);CHKERRQ(ierr);
28845f43ab7SHong Zhang   ierr = PetscBTDestroy(otable);CHKERRQ(ierr);
289632d0f97SHong Zhang 
29045f43ab7SHong Zhang   /* 3. Do local work on this processor's is[] */
29145f43ab7SHong Zhang   /*-------------------------------------------*/
29245f43ab7SHong Zhang   /* make sure there is enough unused space in odata2(=data) array */
29345f43ab7SHong Zhang   len_max = is_max*(Mbs+1); /* max space storing all is[] for this processor */
29410eea884SHong Zhang   if (len_unused < len_max){ /* allocate more space for odata2 */
295240e5896SHong Zhang     ierr = PetscMalloc((len_est+1)*sizeof(int),&odata2);CHKERRQ(ierr);
296240e5896SHong Zhang     odata2_ptr[++nodata2] = odata2;
29710eea884SHong Zhang     len_unused = len_est;
298240e5896SHong Zhang     /* ierr = PetscPrintf(PETSC_COMM_SELF, " [%d] 3. Malloc data2, nodata2: %d\n",rank,nodata2); */
29910eea884SHong Zhang   }
300bfc6803cSHong Zhang 
301240e5896SHong Zhang   data = odata2;
30245f43ab7SHong Zhang   ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,data,table);CHKERRQ(ierr);
30345f43ab7SHong Zhang   ierr = PetscFree(data1_start);CHKERRQ(ierr);
30445f43ab7SHong Zhang 
30545f43ab7SHong Zhang   /* 4. Receive work done on other processors, then merge */
30645f43ab7SHong Zhang   /*------------------------------------------------------*/
30745f43ab7SHong Zhang   /* get max number of messages that this processor expects to recv */
30845f43ab7SHong Zhang   ierr = MPI_Allreduce(len_s,iwork,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
30945f43ab7SHong Zhang   /* ierr = PetscPrintf(PETSC_COMM_SELF," [%d] expects max_len=%d of data2 from others \n",rank,iwork[rank]); */
31045f43ab7SHong Zhang   ierr = PetscMalloc((iwork[rank]+1)*sizeof(int),&data2);CHKERRQ(ierr);
31145f43ab7SHong Zhang   ierr = PetscFree(len_s);CHKERRQ(ierr);
31245f43ab7SHong Zhang 
313632d0f97SHong Zhang   k = 0;
3145483b11dSHong Zhang   while (k < nrqs){
315632d0f97SHong Zhang     /* Receive messages */
316c910923dSHong Zhang     ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status);
317632d0f97SHong Zhang     if (flag){
318bfc6803cSHong Zhang       ierr = MPI_Get_count(&r_status,MPI_INT,&len);CHKERRQ(ierr);
319632d0f97SHong Zhang       proc_id = r_status.MPI_SOURCE;
3200472cc68SHong Zhang       ierr = MPI_Irecv(data2,len,MPI_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr);
3215483b11dSHong Zhang       ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr);
3225483b11dSHong Zhang       /* ierr = PetscPrintf(PETSC_COMM_SELF," [%d] recv %d from [%d], data2:\n",rank,len,proc_id); */
3235483b11dSHong Zhang       if (len > 1+is_max){ /* Add data2 into data */
3240472cc68SHong Zhang         data2_i = data2 + 1 + is_max;
325fc70d14eSHong Zhang         for (i=0; i<is_max; i++){
326fc70d14eSHong Zhang           table_i = table[i];
327bfc6803cSHong Zhang           data_i  = data + 1 + is_max + Mbs*i;
328bfc6803cSHong Zhang           isz     = data[1+i];
3290472cc68SHong Zhang           for (j=0; j<data2[1+i]; j++){
3300472cc68SHong Zhang             col = data2_i[j];
331bfc6803cSHong Zhang             if (!PetscBTLookupSet(table_i,col)) {data_i[isz++] = col;}
332fc70d14eSHong Zhang           }
333bfc6803cSHong Zhang           data[1+i] = isz;
3340472cc68SHong Zhang           if (i < is_max - 1) data2_i += data2[1+i];
335fc70d14eSHong Zhang         }
3365483b11dSHong Zhang       }
337632d0f97SHong Zhang       k++;
338632d0f97SHong Zhang     }
3395483b11dSHong Zhang   }
34045f43ab7SHong Zhang   ierr = PetscFree(data2);CHKERRQ(ierr);
34145f43ab7SHong Zhang   ierr = PetscFree(table);CHKERRQ(ierr);
3425483b11dSHong Zhang 
3435483b11dSHong Zhang   /* phase 1 sends are complete */
3445483b11dSHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&s_status);CHKERRQ(ierr);
3455483b11dSHong Zhang   if (nrqs){
3465483b11dSHong Zhang     ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);
3475483b11dSHong Zhang   }
3485483b11dSHong Zhang   ierr = PetscFree(data1);CHKERRQ(ierr);
3495483b11dSHong Zhang 
350240e5896SHong Zhang   /* phase 2 sends are complete */
3515483b11dSHong Zhang   if (nrqr){
3525483b11dSHong Zhang     ierr = MPI_Waitall(nrqr,s_waits2,s_status);CHKERRQ(ierr);
3535483b11dSHong Zhang   }
35445f43ab7SHong Zhang   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
35545f43ab7SHong Zhang   ierr = PetscFree(s_status);CHKERRQ(ierr);
356632d0f97SHong Zhang 
357c910923dSHong Zhang   /* 5. Create new is[] */
358c910923dSHong Zhang   /*--------------------*/
359c910923dSHong Zhang   for (i=0; i<is_max; i++) {
360bfc6803cSHong Zhang     data_i = data + 1 + is_max + Mbs*i;
361bfc6803cSHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,is+i);CHKERRQ(ierr);
362632d0f97SHong Zhang   }
36345f43ab7SHong Zhang   for (k=0; k<=nodata2; k++){
36445f43ab7SHong Zhang     ierr = PetscFree(odata2_ptr[k]);CHKERRQ(ierr);
36545f43ab7SHong Zhang   }
36645f43ab7SHong Zhang   ierr = PetscFree(odata2_ptr);CHKERRQ(ierr);
3675483b11dSHong Zhang 
368632d0f97SHong Zhang   PetscFunctionReturn(0);
369632d0f97SHong Zhang }
370632d0f97SHong Zhang 
371632d0f97SHong Zhang #undef __FUNCT__
372632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Local"
373632d0f97SHong Zhang /*
374dc008846SHong Zhang    MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do
375632d0f97SHong Zhang        the work on the local processor.
376632d0f97SHong Zhang 
377632d0f97SHong Zhang      Inputs:
378632d0f97SHong Zhang       C      - MAT_MPISBAIJ;
379bfc6803cSHong Zhang       data   - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format.
380bfc6803cSHong Zhang       whose  - whose is[] to be processed,
381bfc6803cSHong Zhang                MINE:  this processor's is[]
382bfc6803cSHong Zhang                OTHER: other processor's is[]
383632d0f97SHong Zhang      Output:
38410eea884SHong Zhang        nidx  - whose = MINE:
3850472cc68SHong Zhang                      holds input and newly found indices in the same format as data
3860472cc68SHong Zhang                whose = OTHER:
3870472cc68SHong Zhang                      only holds the newly found indices
3880472cc68SHong Zhang        table - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'.
389632d0f97SHong Zhang */
39076f244e2SHong Zhang /* Would computation be reduced by swapping the loop 'for each is' and 'for each row'? */
391*6849ba73SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat C,int *data,int whose,int *nidx,PetscBT *table)
392632d0f97SHong Zhang {
393632d0f97SHong Zhang   Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data;
394dc008846SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(c->A)->data;
395dc008846SHong Zhang   Mat_SeqBAIJ  *b = (Mat_SeqBAIJ*)(c->B)->data;
396dfbe8321SBarry Smith   PetscErrorCode ierr;
397dfbe8321SBarry Smith   int row,mbs,Mbs,*nidx_i,col,col_max,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l;
398dc008846SHong Zhang   int          a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n;
399bfc6803cSHong Zhang   PetscBT      table0;  /* mark the indices of input is[] for look up */
400bfc6803cSHong Zhang   PetscBT      table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */
401632d0f97SHong Zhang 
402632d0f97SHong Zhang   PetscFunctionBegin;
40331f99336SHong Zhang   Mbs = c->Mbs; mbs = a->mbs;
404dc008846SHong Zhang   ai = a->i; aj = a->j;
405dc008846SHong Zhang   bi = b->i; bj = b->j;
406632d0f97SHong Zhang   garray = c->garray;
407dc008846SHong Zhang   rstart = c->rstart;
408dc008846SHong Zhang   is_max = data[0];
409632d0f97SHong Zhang 
410fc70d14eSHong Zhang   ierr = PetscBTCreate(Mbs,table0);CHKERRQ(ierr);
411fc70d14eSHong Zhang 
412fc70d14eSHong Zhang   nidx[0] = is_max;
413fc70d14eSHong Zhang   idx_i   = data + is_max + 1; /* ptr to input is[0] array */
414bfc6803cSHong Zhang   nidx_i  = nidx + is_max + 1; /* ptr to output is[0] array */
415dc008846SHong Zhang   for (i=0; i<is_max; i++) { /* for each is */
416dc008846SHong Zhang     isz  = 0;
417fc70d14eSHong Zhang     n = data[1+i]; /* size of input is[i] */
418dc008846SHong Zhang 
41976f244e2SHong Zhang     /* initialize and set table_i(mark idx and nidx) and table0(only mark idx) */
420bfc6803cSHong Zhang     if (whose == MINE){ /* process this processor's is[] */
421bfc6803cSHong Zhang       table_i = table[i];
4220472cc68SHong Zhang       nidx_i  = nidx + 1+ is_max + Mbs*i;
423bfc6803cSHong Zhang     } else {            /* process other processor's is[] - only use one temp table */
424430a0127SHong Zhang       table_i = table[0];
425bfc6803cSHong Zhang     }
426bfc6803cSHong Zhang     ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr);
427bfc6803cSHong Zhang     ierr = PetscBTMemzero(Mbs,table0);CHKERRQ(ierr);
42876f244e2SHong Zhang     if (n==0) {
42976f244e2SHong Zhang        nidx[1+i] = 0; /* size of new is[i] */
43076f244e2SHong Zhang        continue;
43176f244e2SHong Zhang     }
43276f244e2SHong Zhang 
43376f244e2SHong Zhang     isz0 = 0; col_max = 0;
434dc008846SHong Zhang     for (j=0; j<n; j++){
435dc008846SHong Zhang       col = idx_i[j];
436dc008846SHong Zhang       if (col >= Mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"index col %d >= Mbs %d",col,Mbs);
437bfc6803cSHong Zhang       if(!PetscBTLookupSet(table_i,col)) {
438bfc6803cSHong Zhang         ierr = PetscBTSet(table0,col);CHKERRQ(ierr);
4390472cc68SHong Zhang         if (whose == MINE) {nidx_i[isz0] = col;}
440dbe03f88SHong Zhang         if (col_max < col) col_max = col;
441bfc6803cSHong Zhang         isz0++;
442bfc6803cSHong Zhang       }
443632d0f97SHong Zhang     }
444dc008846SHong Zhang 
4450472cc68SHong Zhang     if (whose == MINE) {isz = isz0;}
446fc70d14eSHong Zhang     k = 0;  /* no. of indices from input is[i] that have been examined */
447dc008846SHong Zhang     for (row=0; row<mbs; row++){
448dc008846SHong Zhang       a_start = ai[row]; a_end = ai[row+1];
449dc008846SHong Zhang       b_start = bi[row]; b_end = bi[row+1];
4500472cc68SHong Zhang       if (PetscBTLookup(table0,row+rstart)){ /* row is on input is[i]:
4510472cc68SHong Zhang                                                 do row search: collect all col in this row */
452dc008846SHong Zhang         for (l = a_start; l<a_end ; l++){ /* Amat */
453dc008846SHong Zhang           col = aj[l] + rstart;
454fc70d14eSHong Zhang           if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;}
455dc008846SHong Zhang         }
456dc008846SHong Zhang         for (l = b_start; l<b_end ; l++){ /* Bmat */
457dc008846SHong Zhang           col = garray[bj[l]];
458fc70d14eSHong Zhang           if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;}
459dc008846SHong Zhang         }
460dc008846SHong Zhang         k++;
461dc008846SHong Zhang         if (k >= isz0) break; /* for (row=0; row<mbs; row++) */
4620472cc68SHong Zhang       } else { /* row is not on input is[i]:
4630472cc68SHong Zhang                   do col serach: add row onto nidx_i if there is a col in nidx_i */
464dc008846SHong Zhang         for (l = a_start; l<a_end ; l++){ /* Amat */
465dc008846SHong Zhang           col = aj[l] + rstart;
46676f244e2SHong Zhang           if (col > col_max) break;
467dc008846SHong Zhang           if (PetscBTLookup(table0,col)){
468fc70d14eSHong Zhang             if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;}
469dc008846SHong Zhang             break; /* for l = start; l<end ; l++) */
470632d0f97SHong Zhang           }
471632d0f97SHong Zhang         }
472dc008846SHong Zhang         for (l = b_start; l<b_end ; l++){ /* Bmat */
473dc008846SHong Zhang           col = garray[bj[l]];
47476f244e2SHong Zhang           if (col > col_max) break;
475dc008846SHong Zhang           if (PetscBTLookup(table0,col)){
476fc70d14eSHong Zhang             if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;}
477dc008846SHong Zhang             break; /* for l = start; l<end ; l++) */
478632d0f97SHong Zhang           }
479dc008846SHong Zhang         }
480dc008846SHong Zhang       }
481bfc6803cSHong Zhang     }
482dc008846SHong Zhang 
483dc008846SHong Zhang     if (i < is_max - 1){
484fc70d14eSHong Zhang       idx_i  += n;   /* ptr to input is[i+1] array */
485bfc6803cSHong Zhang       nidx_i += isz; /* ptr to output is[i+1] array */
486dc008846SHong Zhang     }
487fc70d14eSHong Zhang     nidx[1+i] = isz; /* size of new is[i] */
4881ab97a2aSSatish Balay   } /* for each is */
489dc008846SHong Zhang   ierr = PetscBTDestroy(table0);CHKERRQ(ierr);
490632d0f97SHong Zhang 
491632d0f97SHong Zhang   PetscFunctionReturn(0);
492632d0f97SHong Zhang }
493632d0f97SHong Zhang 
494632d0f97SHong Zhang 
495