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 */ 6c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 7c6db04a5SJed Brown #include <petscbt.h> 8632d0f97SHong Zhang 913f74950SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat,PetscInt,IS*); 1013f74950SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat,PetscInt*,PetscInt,PetscInt*,PetscBT*); 11632d0f97SHong Zhang 1213f74950SBarry Smith PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat C,PetscInt is_max,IS is[],PetscInt ov) 13632d0f97SHong Zhang { 146849ba73SBarry Smith PetscErrorCode ierr; 15b11de2a9SHong Zhang PetscInt i,N=C->cmap->N, bs=C->rmap->bs,M=C->rmap->N,Mbs=M/bs,*nidx,isz,iov; 16b11de2a9SHong Zhang IS *is_new,*is_row; 17b11de2a9SHong Zhang Mat *submats; 18b11de2a9SHong Zhang Mat_MPISBAIJ *c=(Mat_MPISBAIJ*)C->data; 19b11de2a9SHong Zhang Mat_SeqSBAIJ *asub_i; 20b11de2a9SHong Zhang PetscBT table; 21b11de2a9SHong Zhang PetscInt *ai,brow,nz,nis,l,nmax,nstages_local,nstages,max_no,pos; 22b11de2a9SHong Zhang const PetscInt *idx; 23b829d874SHong Zhang PetscBool flg; 24632d0f97SHong Zhang 25632d0f97SHong Zhang PetscFunctionBegin; 26785e854fSJed Brown ierr = PetscMalloc1(is_max,&is_new);CHKERRQ(ierr); 27632d0f97SHong Zhang /* Convert the indices into block format */ 2805d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,is_max,is,is_new);CHKERRQ(ierr); 2926fbe8dcSKarl Rupp if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n"); 30db41cccfSHong Zhang 31b11de2a9SHong Zhang /* ----- previous non-scalable implementation ----- */ 32a0769a71SSatish Balay flg = PETSC_FALSE; 33c5929fdfSBarry Smith ierr = PetscOptionsHasName(NULL,NULL, "-IncreaseOverlap_old", &flg);CHKERRQ(ierr); 347a868f3eSHong Zhang if (flg) { /* previous non-scalable implementation */ 357a868f3eSHong Zhang printf("use previous non-scalable implementation...\n"); 36632d0f97SHong Zhang for (i=0; i<ov; ++i) { 37c910923dSHong Zhang ierr = MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr); 38632d0f97SHong Zhang } 391c4982ebSHong Zhang } else { /* implementation using modified BAIJ routines */ 401c4982ebSHong Zhang 41854ce69bSBarry Smith ierr = PetscMalloc1(Mbs+1,&nidx);CHKERRQ(ierr); 4253b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&table);CHKERRQ(ierr); /* for column search */ 43db41cccfSHong Zhang 44db41cccfSHong Zhang /* Create is_row */ 45785e854fSJed Brown ierr = PetscMalloc1(is_max,&is_row);CHKERRQ(ierr); 464da72fa9SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,Mbs,0,1,&is_row[0]);CHKERRQ(ierr); 4726fbe8dcSKarl Rupp 484da72fa9SHong Zhang for (i=1; i<is_max; i++) { 494da72fa9SHong Zhang is_row[i] = is_row[0]; /* reuse is_row[0] */ 504da72fa9SHong Zhang } 51db41cccfSHong Zhang 527dae84e0SHong Zhang /* Allocate memory to hold all the submatrices - Modified from MatCreateSubMatrices_MPIBAIJ() */ 53854ce69bSBarry Smith ierr = PetscMalloc1(is_max+1,&submats);CHKERRQ(ierr); 54b11de2a9SHong Zhang 55db41cccfSHong Zhang /* Determine the number of stages through which submatrices are done */ 56db41cccfSHong Zhang nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 57db41cccfSHong Zhang if (!nmax) nmax = 1; 58db41cccfSHong Zhang nstages_local = is_max/nmax + ((is_max % nmax) ? 1 : 0); 59db41cccfSHong Zhang 60db41cccfSHong Zhang /* Make sure every processor loops through the nstages */ 61820f2d46SBarry Smith ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRMPI(ierr); 62b11de2a9SHong Zhang 63b11de2a9SHong Zhang for (iov=0; iov<ov; ++iov) { 64b11de2a9SHong Zhang /* 1) Get submats for column search */ 65db41cccfSHong Zhang for (i=0,pos=0; i<nstages; i++) { 66db41cccfSHong Zhang if (pos+nmax <= is_max) max_no = nmax; 67db41cccfSHong Zhang else if (pos == is_max) max_no = 0; 68db41cccfSHong Zhang else max_no = is_max-pos; 69b829d874SHong Zhang c->ijonly = PETSC_TRUE; /* only matrix data structures are requested */ 70476417e5SBarry Smith /* The resulting submatrices should be BAIJ, not SBAIJ, hence we change this value to trigger that */ 71476417e5SBarry Smith ierr = PetscStrcpy(((PetscObject)c->A)->type_name,MATSEQBAIJ);CHKERRQ(ierr); 727dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIBAIJ_local(C,max_no,is_row+pos,is_new+pos,MAT_INITIAL_MATRIX,submats+pos);CHKERRQ(ierr); 73476417e5SBarry Smith ierr = PetscStrcpy(((PetscObject)c->A)->type_name,MATSEQSBAIJ);CHKERRQ(ierr); 74db41cccfSHong Zhang pos += max_no; 75db41cccfSHong Zhang } 76db41cccfSHong Zhang 77b11de2a9SHong Zhang /* 2) Row search */ 78db41cccfSHong Zhang ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr); 79db41cccfSHong Zhang 80b11de2a9SHong Zhang /* 3) Column search */ 81db41cccfSHong Zhang for (i=0; i<is_max; i++) { 82db41cccfSHong Zhang asub_i = (Mat_SeqSBAIJ*)submats[i]->data; 83feb237baSPierre Jolivet ai = asub_i->i; 84db41cccfSHong Zhang 85db41cccfSHong Zhang /* put is_new obtained from MatIncreaseOverlap_MPIBAIJ() to table */ 86db41cccfSHong Zhang ierr = PetscBTMemzero(Mbs,table);CHKERRQ(ierr); 87db41cccfSHong Zhang 88db41cccfSHong Zhang ierr = ISGetIndices(is_new[i],&idx);CHKERRQ(ierr); 89db41cccfSHong Zhang ierr = ISGetLocalSize(is_new[i],&nis);CHKERRQ(ierr); 90db41cccfSHong Zhang for (l=0; l<nis; l++) { 91db41cccfSHong Zhang ierr = PetscBTSet(table,idx[l]);CHKERRQ(ierr); 92db41cccfSHong Zhang nidx[l] = idx[l]; 93db41cccfSHong Zhang } 94db41cccfSHong Zhang isz = nis; 95db41cccfSHong Zhang 96b11de2a9SHong Zhang /* add column entries to table */ 97db41cccfSHong Zhang for (brow=0; brow<Mbs; brow++) { 98db41cccfSHong Zhang nz = ai[brow+1] - ai[brow]; 99db41cccfSHong Zhang if (nz) { 100db41cccfSHong Zhang if (!PetscBTLookupSet(table,brow)) nidx[isz++] = brow; 101db41cccfSHong Zhang } 102db41cccfSHong Zhang } 103db41cccfSHong Zhang ierr = ISRestoreIndices(is_new[i],&idx);CHKERRQ(ierr); 1046bf464f9SBarry Smith ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr); 105db41cccfSHong Zhang 106db41cccfSHong Zhang /* create updated is_new */ 107db41cccfSHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,is_new+i);CHKERRQ(ierr); 108db41cccfSHong Zhang } 109db41cccfSHong Zhang 110db41cccfSHong Zhang /* Free tmp spaces */ 111db41cccfSHong Zhang for (i=0; i<is_max; i++) { 1126bf464f9SBarry Smith ierr = MatDestroy(&submats[i]);CHKERRQ(ierr); 113db41cccfSHong Zhang } 114db41cccfSHong Zhang } 115b829d874SHong Zhang 11694bacf5dSBarry Smith ierr = PetscBTDestroy(&table);CHKERRQ(ierr); 117db41cccfSHong Zhang ierr = PetscFree(submats);CHKERRQ(ierr); 1186bf464f9SBarry Smith ierr = ISDestroy(&is_row[0]);CHKERRQ(ierr); 119db41cccfSHong Zhang ierr = PetscFree(is_row);CHKERRQ(ierr); 120db41cccfSHong Zhang ierr = PetscFree(nidx);CHKERRQ(ierr); 121b11de2a9SHong Zhang 122db41cccfSHong Zhang } 123db41cccfSHong Zhang 1246bf464f9SBarry Smith for (i=0; i<is_max; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 12505d8c843SHong Zhang ierr = ISExpandIndicesGeneral(N,N,bs,is_max,is_new,is);CHKERRQ(ierr); 126db41cccfSHong Zhang 1276bf464f9SBarry Smith for (i=0; i<is_max; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 128632d0f97SHong Zhang ierr = PetscFree(is_new);CHKERRQ(ierr); 129632d0f97SHong Zhang PetscFunctionReturn(0); 130632d0f97SHong Zhang } 131632d0f97SHong Zhang 1324a69c536SBarry Smith typedef enum {MINE,OTHER} WhoseOwner; 1330472cc68SHong Zhang /* data1, odata1 and odata2 are packed in the format (for communication): 134a2a9f22aSHong Zhang data[0] = is_max, no of is 135a2a9f22aSHong Zhang data[1] = size of is[0] 136a2a9f22aSHong Zhang ... 137a2a9f22aSHong Zhang data[is_max] = size of is[is_max-1] 138a2a9f22aSHong Zhang data[is_max + 1] = data(is[0]) 139a2a9f22aSHong Zhang ... 140a2a9f22aSHong Zhang data[is_max+1+sum(size of is[k]), k=0,...,i-1] = data(is[i]) 141a2a9f22aSHong Zhang ... 1420472cc68SHong Zhang data2 is packed in the format (for creating output is[]): 1430472cc68SHong Zhang data[0] = is_max, no of is 1440472cc68SHong Zhang data[1] = size of is[0] 1450472cc68SHong Zhang ... 1460472cc68SHong Zhang data[is_max] = size of is[is_max-1] 1470472cc68SHong Zhang data[is_max + 1] = data(is[0]) 1480472cc68SHong Zhang ... 1490472cc68SHong Zhang data[is_max + 1 + Mbs*i) = data(is[i]) 1500472cc68SHong Zhang ... 151a2a9f22aSHong Zhang */ 15213f74950SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat C,PetscInt is_max,IS is[]) 153632d0f97SHong Zhang { 154632d0f97SHong Zhang Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 1556849ba73SBarry Smith PetscErrorCode ierr; 156adcec1e5SJed Brown PetscMPIInt size,rank,tag1,tag2,*len_s,nrqr,nrqs,*id_r1,*len_r1,flag,len,*iwork; 1575d0c19d7SBarry Smith const PetscInt *idx_i; 15826fbe8dcSKarl Rupp PetscInt idx,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i; 15926fbe8dcSKarl Rupp PetscInt Mbs,i,j,k,*odata1,*odata2; 160f4259b30SLisandro Dalcin PetscInt proc_id,**odata2_ptr,*ctable=NULL,*btable,len_max,len_est; 161adcec1e5SJed Brown PetscInt proc_end=0,len_unused,nodata2; 16213f74950SBarry Smith PetscInt ois_max; /* max no of is[] in each of processor */ 163bfc6803cSHong Zhang char *t_p; 164632d0f97SHong Zhang MPI_Comm comm; 165e8527bf2SHong Zhang MPI_Request *s_waits1,*s_waits2,r_req; 166632d0f97SHong Zhang MPI_Status *s_status,r_status; 16745f43ab7SHong Zhang PetscBT *table; /* mark indices of this processor's is[] */ 168fc70d14eSHong Zhang PetscBT table_i; 169bfc6803cSHong Zhang PetscBT otable; /* mark indices of other processors' is[] */ 170d0f46423SBarry Smith PetscInt bs=C->rmap->bs,Bn = c->B->cmap->n,Bnbs = Bn/bs,*Bowners; 17110eea884SHong Zhang IS garray_local,garray_gl; 1725483b11dSHong Zhang 173632d0f97SHong Zhang PetscFunctionBegin; 174ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 175632d0f97SHong Zhang size = c->size; 176632d0f97SHong Zhang rank = c->rank; 177632d0f97SHong Zhang Mbs = c->Mbs; 178632d0f97SHong Zhang 179c910923dSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 180c910923dSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 181c910923dSHong Zhang 182430a0127SHong Zhang /* create tables used in 183430a0127SHong Zhang step 1: table[i] - mark c->garray of proc [i] 18445f43ab7SHong Zhang step 3: table[i] - mark indices of is[i] when whose=MINE 185430a0127SHong Zhang table[0] - mark incideces of is[] when whose=OTHER */ 186*2f613bf5SBarry Smith len = PetscMax(is_max, size); 187dcca6d9dSJed Brown ierr = PetscMalloc2(len,&table,(Mbs/PETSC_BITS_PER_BYTE+1)*len,&t_p);CHKERRQ(ierr); 188430a0127SHong Zhang for (i=0; i<len; i++) { 189430a0127SHong Zhang table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 190430a0127SHong Zhang } 191430a0127SHong Zhang 192820f2d46SBarry Smith ierr = MPIU_Allreduce(&is_max,&ois_max,1,MPIU_INT,MPI_MAX,comm);CHKERRMPI(ierr); 19310eea884SHong Zhang 1945483b11dSHong Zhang /* 1. Send this processor's is[] to other processors */ 1955483b11dSHong Zhang /*---------------------------------------------------*/ 196e8527bf2SHong Zhang /* allocate spaces */ 197785e854fSJed Brown ierr = PetscMalloc1(is_max,&n);CHKERRQ(ierr); 1985483b11dSHong Zhang len = 0; 199c910923dSHong Zhang for (i=0; i<is_max; i++) { 200632d0f97SHong Zhang ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 201632d0f97SHong Zhang len += n[i]; 202632d0f97SHong Zhang } 203958c9bccSBarry Smith if (!len) { 2045483b11dSHong Zhang is_max = 0; 2055483b11dSHong Zhang } else { 2065483b11dSHong Zhang len += 1 + is_max; /* max length of data1 for one processor */ 2075483b11dSHong Zhang } 2085483b11dSHong Zhang 209854ce69bSBarry Smith ierr = PetscMalloc1(size*len+1,&data1);CHKERRQ(ierr); 210785e854fSJed Brown ierr = PetscMalloc1(size,&data1_start);CHKERRQ(ierr); 2115483b11dSHong Zhang for (i=0; i<size; i++) data1_start[i] = data1 + i*len; 2125483b11dSHong Zhang 213dcca6d9dSJed Brown ierr = PetscMalloc4(size,&len_s,size,&btable,size,&iwork,size+1,&Bowners);CHKERRQ(ierr); 214e8527bf2SHong Zhang 21576f244e2SHong Zhang /* gather c->garray from all processors */ 21670b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,Bnbs,c->garray,PETSC_COPY_VALUES,&garray_local);CHKERRQ(ierr); 21776f244e2SHong Zhang ierr = ISAllGather(garray_local, &garray_gl);CHKERRQ(ierr); 2186bf464f9SBarry Smith ierr = ISDestroy(&garray_local);CHKERRQ(ierr); 219ffc4695bSBarry Smith ierr = MPI_Allgather(&Bnbs,1,MPIU_INT,Bowners+1,1,MPIU_INT,comm);CHKERRMPI(ierr); 22026fbe8dcSKarl Rupp 22176f244e2SHong Zhang Bowners[0] = 0; 22276f244e2SHong Zhang for (i=0; i<size; i++) Bowners[i+1] += Bowners[i]; 22376f244e2SHong Zhang 2245483b11dSHong Zhang if (is_max) { 225430a0127SHong Zhang /* hash table ctable which maps c->row to proc_id) */ 226785e854fSJed Brown ierr = PetscMalloc1(Mbs,&ctable);CHKERRQ(ierr); 2275483b11dSHong Zhang for (proc_id=0,j=0; proc_id<size; proc_id++) { 22826fbe8dcSKarl Rupp for (; j<C->rmap->range[proc_id+1]/bs; j++) ctable[j] = proc_id; 2295483b11dSHong Zhang } 2305483b11dSHong Zhang 231430a0127SHong Zhang /* hash tables marking c->garray */ 23222d28d08SBarry Smith ierr = ISGetIndices(garray_gl,&idx_i);CHKERRQ(ierr); 2335483b11dSHong Zhang for (i=0; i<size; i++) { 234430a0127SHong Zhang table_i = table[i]; 235430a0127SHong Zhang ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); 236430a0127SHong Zhang for (j = Bowners[i]; j<Bowners[i+1]; j++) { /* go through B cols of proc[i]*/ 237430a0127SHong Zhang ierr = PetscBTSet(table_i,idx_i[j]);CHKERRQ(ierr); 2385483b11dSHong Zhang } 2395483b11dSHong Zhang } 24010eea884SHong Zhang ierr = ISRestoreIndices(garray_gl,&idx_i);CHKERRQ(ierr); 2415483b11dSHong Zhang } /* if (is_max) */ 2426bf464f9SBarry Smith ierr = ISDestroy(&garray_gl);CHKERRQ(ierr); 2435483b11dSHong Zhang 2445483b11dSHong Zhang /* evaluate communication - mesg to who, length, and buffer space */ 245e8527bf2SHong Zhang for (i=0; i<size; i++) len_s[i] = 0; 2465483b11dSHong Zhang 2475483b11dSHong Zhang /* header of data1 */ 2485483b11dSHong Zhang for (proc_id=0; proc_id<size; proc_id++) { 2495483b11dSHong Zhang iwork[proc_id] = 0; 2505483b11dSHong Zhang *data1_start[proc_id] = is_max; 2515483b11dSHong Zhang data1_start[proc_id]++; 2525483b11dSHong Zhang for (j=0; j<is_max; j++) { 2535483b11dSHong Zhang if (proc_id == rank) { 2545483b11dSHong Zhang *data1_start[proc_id] = n[j]; 2555483b11dSHong Zhang } else { 2565483b11dSHong Zhang *data1_start[proc_id] = 0; 2575483b11dSHong Zhang } 2585483b11dSHong Zhang data1_start[proc_id]++; 2595483b11dSHong Zhang } 2605483b11dSHong Zhang } 2615483b11dSHong Zhang 2625483b11dSHong Zhang for (i=0; i<is_max; i++) { 2635483b11dSHong Zhang ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr); 2645483b11dSHong Zhang for (j=0; j<n[i]; j++) { 2655483b11dSHong Zhang idx = idx_i[j]; 2665483b11dSHong Zhang *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */ 2675483b11dSHong Zhang proc_end = ctable[idx]; 2685483b11dSHong Zhang for (proc_id=0; proc_id<=proc_end; proc_id++) { /* for others to process */ 269e8527bf2SHong Zhang if (proc_id == rank) continue; /* done before this loop */ 27026fbe8dcSKarl Rupp if (proc_id < proc_end && !PetscBTLookup(table[proc_id],idx)) continue; /* no need for sending idx to [proc_id] */ 2715483b11dSHong Zhang *data1_start[proc_id] = idx; data1_start[proc_id]++; 2725483b11dSHong Zhang len_s[proc_id]++; 2735483b11dSHong Zhang } 2745483b11dSHong Zhang } 2755483b11dSHong Zhang /* update header data */ 2762cfbe0a4SHong Zhang for (proc_id=0; proc_id<size; proc_id++) { 2775483b11dSHong Zhang if (proc_id== rank) continue; 2785483b11dSHong Zhang *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id]; 2795483b11dSHong Zhang iwork[proc_id] = len_s[proc_id]; 2805483b11dSHong Zhang } 2815483b11dSHong Zhang ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr); 282e8527bf2SHong Zhang } 2835483b11dSHong Zhang 2845483b11dSHong Zhang nrqs = 0; nrqr = 0; 2855483b11dSHong Zhang for (i=0; i<size; i++) { 2865483b11dSHong Zhang data1_start[i] = data1 + i*len; 2875483b11dSHong Zhang if (len_s[i]) { 2885483b11dSHong Zhang nrqs++; 2895483b11dSHong Zhang len_s[i] += 1 + is_max; /* add no. of header msg */ 2905483b11dSHong Zhang } 2915483b11dSHong Zhang } 2925483b11dSHong Zhang 2935483b11dSHong Zhang for (i=0; i<is_max; i++) { 2946bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 295632d0f97SHong Zhang } 296bfc6803cSHong Zhang ierr = PetscFree(n);CHKERRQ(ierr); 29705b42c5fSBarry Smith ierr = PetscFree(ctable);CHKERRQ(ierr); 2985483b11dSHong Zhang 2995483b11dSHong Zhang /* Determine the number of messages to expect, their lengths, from from-ids */ 3000298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrqr);CHKERRQ(ierr); 3015483b11dSHong Zhang ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);CHKERRQ(ierr); 302632d0f97SHong Zhang 303632d0f97SHong Zhang /* Now post the sends */ 304dcca6d9dSJed Brown ierr = PetscMalloc2(size,&s_waits1,size,&s_waits2);CHKERRQ(ierr); 305632d0f97SHong Zhang k = 0; 3065483b11dSHong Zhang for (proc_id=0; proc_id<size; proc_id++) { /* send data1 to processor [proc_id] */ 3075483b11dSHong Zhang if (len_s[proc_id]) { 308ffc4695bSBarry Smith ierr = MPI_Isend(data1_start[proc_id],len_s[proc_id],MPIU_INT,proc_id,tag1,comm,s_waits1+k);CHKERRMPI(ierr); 309632d0f97SHong Zhang k++; 310632d0f97SHong Zhang } 311632d0f97SHong Zhang } 312632d0f97SHong Zhang 31345f43ab7SHong Zhang /* 2. Receive other's is[] and process. Then send back */ 314bfc6803cSHong Zhang /*-----------------------------------------------------*/ 3155483b11dSHong Zhang len = 0; 3165483b11dSHong Zhang for (i=0; i<nrqr; i++) { 3175483b11dSHong Zhang if (len_r1[i] > len) len = len_r1[i]; 3185483b11dSHong Zhang } 31945f43ab7SHong Zhang ierr = PetscFree(len_r1);CHKERRQ(ierr); 32045f43ab7SHong Zhang ierr = PetscFree(id_r1);CHKERRQ(ierr); 32145f43ab7SHong Zhang 32226fbe8dcSKarl Rupp for (proc_id=0; proc_id<size; proc_id++) len_s[proc_id] = iwork[proc_id] = 0; 32345f43ab7SHong Zhang 324854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&odata1);CHKERRQ(ierr); 325785e854fSJed Brown ierr = PetscMalloc1(size,&odata2_ptr);CHKERRQ(ierr); 32653b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&otable);CHKERRQ(ierr); 32710eea884SHong Zhang 32810eea884SHong Zhang len_max = ois_max*(Mbs+1); /* max space storing all is[] for each receive */ 329240e5896SHong Zhang len_est = 2*len_max; /* estimated space of storing is[] for all receiving messages */ 330854ce69bSBarry Smith ierr = PetscMalloc1(len_est+1,&odata2);CHKERRQ(ierr); 33110eea884SHong Zhang nodata2 = 0; /* nodata2+1: num of PetscMalloc(,&odata2_ptr[]) called */ 33226fbe8dcSKarl Rupp 333240e5896SHong Zhang odata2_ptr[nodata2] = odata2; 33426fbe8dcSKarl Rupp 33510eea884SHong Zhang len_unused = len_est; /* unused space in the array odata2_ptr[nodata2]-- needs to be >= len_max */ 33610eea884SHong Zhang 337632d0f97SHong Zhang k = 0; 3385483b11dSHong Zhang while (k < nrqr) { 339632d0f97SHong Zhang /* Receive messages */ 340ffc4695bSBarry Smith ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);CHKERRMPI(ierr); 341632d0f97SHong Zhang if (flag) { 34255b25c41SPierre Jolivet ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRMPI(ierr); 343632d0f97SHong Zhang proc_id = r_status.MPI_SOURCE; 34455b25c41SPierre Jolivet ierr = MPI_Irecv(odata1,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRMPI(ierr); 34555b25c41SPierre Jolivet ierr = MPI_Wait(&r_req,&r_status);CHKERRMPI(ierr); 346632d0f97SHong Zhang 347fc70d14eSHong Zhang /* Process messages */ 348240e5896SHong Zhang /* make sure there is enough unused space in odata2 array */ 34910eea884SHong Zhang if (len_unused < len_max) { /* allocate more space for odata2 */ 350854ce69bSBarry Smith ierr = PetscMalloc1(len_est+1,&odata2);CHKERRQ(ierr); 35126fbe8dcSKarl Rupp 352240e5896SHong Zhang odata2_ptr[++nodata2] = odata2; 35326fbe8dcSKarl Rupp 35410eea884SHong Zhang len_unused = len_est; 35510eea884SHong Zhang } 35610eea884SHong Zhang 35710eea884SHong Zhang ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,odata2,&otable);CHKERRQ(ierr); 358a2a9f22aSHong Zhang len = 1 + odata2[0]; 35926fbe8dcSKarl Rupp for (i=0; i<odata2[0]; i++) len += odata2[1 + i]; 360632d0f97SHong Zhang 361632d0f97SHong Zhang /* Send messages back */ 362ffc4695bSBarry Smith ierr = MPI_Isend(odata2,len,MPIU_INT,proc_id,tag2,comm,s_waits2+k);CHKERRMPI(ierr); 363fc70d14eSHong Zhang k++; 36410eea884SHong Zhang odata2 += len; 36510eea884SHong Zhang len_unused -= len; 36645f43ab7SHong Zhang len_s[proc_id] = len; /* num of messages sending back to [proc_id] by this proc */ 367632d0f97SHong Zhang } 3685483b11dSHong Zhang } 3695483b11dSHong Zhang ierr = PetscFree(odata1);CHKERRQ(ierr); 37094bacf5dSBarry Smith ierr = PetscBTDestroy(&otable);CHKERRQ(ierr); 371632d0f97SHong Zhang 37245f43ab7SHong Zhang /* 3. Do local work on this processor's is[] */ 37345f43ab7SHong Zhang /*-------------------------------------------*/ 37445f43ab7SHong Zhang /* make sure there is enough unused space in odata2(=data) array */ 37545f43ab7SHong Zhang len_max = is_max*(Mbs+1); /* max space storing all is[] for this processor */ 37610eea884SHong Zhang if (len_unused < len_max) { /* allocate more space for odata2 */ 377854ce69bSBarry Smith ierr = PetscMalloc1(len_est+1,&odata2);CHKERRQ(ierr); 37826fbe8dcSKarl Rupp 379240e5896SHong Zhang odata2_ptr[++nodata2] = odata2; 38010eea884SHong Zhang } 381bfc6803cSHong Zhang 382240e5896SHong Zhang data = odata2; 38345f43ab7SHong Zhang ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,data,table);CHKERRQ(ierr); 38445f43ab7SHong Zhang ierr = PetscFree(data1_start);CHKERRQ(ierr); 38545f43ab7SHong Zhang 38645f43ab7SHong Zhang /* 4. Receive work done on other processors, then merge */ 38745f43ab7SHong Zhang /*------------------------------------------------------*/ 38845f43ab7SHong Zhang /* get max number of messages that this processor expects to recv */ 389820f2d46SBarry Smith ierr = MPIU_Allreduce(len_s,iwork,size,MPI_INT,MPI_MAX,comm);CHKERRMPI(ierr); 390854ce69bSBarry Smith ierr = PetscMalloc1(iwork[rank]+1,&data2);CHKERRQ(ierr); 39174ed9c26SBarry Smith ierr = PetscFree4(len_s,btable,iwork,Bowners);CHKERRQ(ierr); 39245f43ab7SHong Zhang 393632d0f97SHong Zhang k = 0; 3945483b11dSHong Zhang while (k < nrqs) { 395632d0f97SHong Zhang /* Receive messages */ 396ffc4695bSBarry Smith ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status);CHKERRMPI(ierr); 397632d0f97SHong Zhang if (flag) { 398ffc4695bSBarry Smith ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRMPI(ierr); 39926fbe8dcSKarl Rupp 400632d0f97SHong Zhang proc_id = r_status.MPI_SOURCE; 40126fbe8dcSKarl Rupp 402ffc4695bSBarry Smith ierr = MPI_Irecv(data2,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRMPI(ierr); 403ffc4695bSBarry Smith ierr = MPI_Wait(&r_req,&r_status);CHKERRMPI(ierr); 4045483b11dSHong Zhang if (len > 1+is_max) { /* Add data2 into data */ 4050472cc68SHong Zhang data2_i = data2 + 1 + is_max; 406fc70d14eSHong Zhang for (i=0; i<is_max; i++) { 407fc70d14eSHong Zhang table_i = table[i]; 408bfc6803cSHong Zhang data_i = data + 1 + is_max + Mbs*i; 409bfc6803cSHong Zhang isz = data[1+i]; 4100472cc68SHong Zhang for (j=0; j<data2[1+i]; j++) { 4110472cc68SHong Zhang col = data2_i[j]; 41226fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,col)) data_i[isz++] = col; 413fc70d14eSHong Zhang } 414bfc6803cSHong Zhang data[1+i] = isz; 4150472cc68SHong Zhang if (i < is_max - 1) data2_i += data2[1+i]; 416fc70d14eSHong Zhang } 4175483b11dSHong Zhang } 418632d0f97SHong Zhang k++; 419632d0f97SHong Zhang } 4205483b11dSHong Zhang } 42145f43ab7SHong Zhang ierr = PetscFree(data2);CHKERRQ(ierr); 42274ed9c26SBarry Smith ierr = PetscFree2(table,t_p);CHKERRQ(ierr); 4235483b11dSHong Zhang 4245483b11dSHong Zhang /* phase 1 sends are complete */ 425785e854fSJed Brown ierr = PetscMalloc1(size,&s_status);CHKERRQ(ierr); 426ffc4695bSBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRMPI(ierr);} 4275483b11dSHong Zhang ierr = PetscFree(data1);CHKERRQ(ierr); 4285483b11dSHong Zhang 429240e5896SHong Zhang /* phase 2 sends are complete */ 430ffc4695bSBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status);CHKERRMPI(ierr);} 43174ed9c26SBarry Smith ierr = PetscFree2(s_waits1,s_waits2);CHKERRQ(ierr); 43245f43ab7SHong Zhang ierr = PetscFree(s_status);CHKERRQ(ierr); 433632d0f97SHong Zhang 434c910923dSHong Zhang /* 5. Create new is[] */ 435c910923dSHong Zhang /*--------------------*/ 436c910923dSHong Zhang for (i=0; i<is_max; i++) { 437bfc6803cSHong Zhang data_i = data + 1 + is_max + Mbs*i; 43870b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 439632d0f97SHong Zhang } 44045f43ab7SHong Zhang for (k=0; k<=nodata2; k++) { 44145f43ab7SHong Zhang ierr = PetscFree(odata2_ptr[k]);CHKERRQ(ierr); 44245f43ab7SHong Zhang } 44345f43ab7SHong Zhang ierr = PetscFree(odata2_ptr);CHKERRQ(ierr); 444632d0f97SHong Zhang PetscFunctionReturn(0); 445632d0f97SHong Zhang } 446632d0f97SHong Zhang 447632d0f97SHong Zhang /* 448dc008846SHong Zhang MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do 449632d0f97SHong Zhang the work on the local processor. 450632d0f97SHong Zhang 451632d0f97SHong Zhang Inputs: 452632d0f97SHong Zhang C - MAT_MPISBAIJ; 453bfc6803cSHong Zhang data - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format. 454bfc6803cSHong Zhang whose - whose is[] to be processed, 455bfc6803cSHong Zhang MINE: this processor's is[] 456bfc6803cSHong Zhang OTHER: other processor's is[] 457632d0f97SHong Zhang Output: 45810eea884SHong Zhang nidx - whose = MINE: 4590472cc68SHong Zhang holds input and newly found indices in the same format as data 4600472cc68SHong Zhang whose = OTHER: 4610472cc68SHong Zhang only holds the newly found indices 4620472cc68SHong Zhang table - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'. 463632d0f97SHong Zhang */ 46476f244e2SHong Zhang /* Would computation be reduced by swapping the loop 'for each is' and 'for each row'? */ 46513f74950SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat C,PetscInt *data,PetscInt whose,PetscInt *nidx,PetscBT *table) 466632d0f97SHong Zhang { 467632d0f97SHong Zhang Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 468dc008846SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(c->A)->data; 469dc008846SHong Zhang Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(c->B)->data; 470dfbe8321SBarry Smith PetscErrorCode ierr; 47113f74950SBarry Smith PetscInt row,mbs,Mbs,*nidx_i,col,col_max,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l; 47213f74950SBarry Smith PetscInt a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n; 473bfc6803cSHong Zhang PetscBT table0; /* mark the indices of input is[] for look up */ 474bfc6803cSHong Zhang PetscBT table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */ 475632d0f97SHong Zhang 476632d0f97SHong Zhang PetscFunctionBegin; 47731f99336SHong Zhang Mbs = c->Mbs; mbs = a->mbs; 478dc008846SHong Zhang ai = a->i; aj = a->j; 479dc008846SHong Zhang bi = b->i; bj = b->j; 480632d0f97SHong Zhang garray = c->garray; 481899cda47SBarry Smith rstart = c->rstartbs; 482dc008846SHong Zhang is_max = data[0]; 483632d0f97SHong Zhang 48453b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&table0);CHKERRQ(ierr); 485fc70d14eSHong Zhang 486fc70d14eSHong Zhang nidx[0] = is_max; 487fc70d14eSHong Zhang idx_i = data + is_max + 1; /* ptr to input is[0] array */ 488bfc6803cSHong Zhang nidx_i = nidx + is_max + 1; /* ptr to output is[0] array */ 489dc008846SHong Zhang for (i=0; i<is_max; i++) { /* for each is */ 490dc008846SHong Zhang isz = 0; 491fc70d14eSHong Zhang n = data[1+i]; /* size of input is[i] */ 492dc008846SHong Zhang 49376f244e2SHong Zhang /* initialize and set table_i(mark idx and nidx) and table0(only mark idx) */ 494bfc6803cSHong Zhang if (whose == MINE) { /* process this processor's is[] */ 495bfc6803cSHong Zhang table_i = table[i]; 4960472cc68SHong Zhang nidx_i = nidx + 1+ is_max + Mbs*i; 497bfc6803cSHong Zhang } else { /* process other processor's is[] - only use one temp table */ 498430a0127SHong Zhang table_i = table[0]; 499bfc6803cSHong Zhang } 500bfc6803cSHong Zhang ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); 501bfc6803cSHong Zhang ierr = PetscBTMemzero(Mbs,table0);CHKERRQ(ierr); 50276f244e2SHong Zhang if (n==0) { 50376f244e2SHong Zhang nidx[1+i] = 0; /* size of new is[i] */ 50476f244e2SHong Zhang continue; 50576f244e2SHong Zhang } 50676f244e2SHong Zhang 50776f244e2SHong Zhang isz0 = 0; col_max = 0; 508dc008846SHong Zhang for (j=0; j<n; j++) { 509dc008846SHong Zhang col = idx_i[j]; 510e32f2f54SBarry Smith if (col >= Mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index col %D >= Mbs %D",col,Mbs); 511bfc6803cSHong Zhang if (!PetscBTLookupSet(table_i,col)) { 512bfc6803cSHong Zhang ierr = PetscBTSet(table0,col);CHKERRQ(ierr); 51326fbe8dcSKarl Rupp if (whose == MINE) nidx_i[isz0] = col; 514dbe03f88SHong Zhang if (col_max < col) col_max = col; 515bfc6803cSHong Zhang isz0++; 516bfc6803cSHong Zhang } 517632d0f97SHong Zhang } 518dc008846SHong Zhang 51926fbe8dcSKarl Rupp if (whose == MINE) isz = isz0; 520fc70d14eSHong Zhang k = 0; /* no. of indices from input is[i] that have been examined */ 521dc008846SHong Zhang for (row=0; row<mbs; row++) { 522dc008846SHong Zhang a_start = ai[row]; a_end = ai[row+1]; 523dc008846SHong Zhang b_start = bi[row]; b_end = bi[row+1]; 5240472cc68SHong Zhang if (PetscBTLookup(table0,row+rstart)) { /* row is on input is[i]: 5250472cc68SHong Zhang do row search: collect all col in this row */ 526dc008846SHong Zhang for (l = a_start; l<a_end ; l++) { /* Amat */ 527dc008846SHong Zhang col = aj[l] + rstart; 52826fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,col)) nidx_i[isz++] = col; 529dc008846SHong Zhang } 530dc008846SHong Zhang for (l = b_start; l<b_end ; l++) { /* Bmat */ 531dc008846SHong Zhang col = garray[bj[l]]; 53226fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,col)) nidx_i[isz++] = col; 533dc008846SHong Zhang } 534dc008846SHong Zhang k++; 535dc008846SHong Zhang if (k >= isz0) break; /* for (row=0; row<mbs; row++) */ 5360472cc68SHong Zhang } else { /* row is not on input is[i]: 5370472cc68SHong Zhang do col serach: add row onto nidx_i if there is a col in nidx_i */ 538dc008846SHong Zhang for (l = a_start; l<a_end; l++) { /* Amat */ 539dc008846SHong Zhang col = aj[l] + rstart; 54076f244e2SHong Zhang if (col > col_max) break; 541dc008846SHong Zhang if (PetscBTLookup(table0,col)) { 54226fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row+rstart)) nidx_i[isz++] = row+rstart; 543dc008846SHong Zhang break; /* for l = start; l<end ; l++) */ 544632d0f97SHong Zhang } 545632d0f97SHong Zhang } 546dc008846SHong Zhang for (l = b_start; l<b_end; l++) { /* Bmat */ 547dc008846SHong Zhang col = garray[bj[l]]; 54876f244e2SHong Zhang if (col > col_max) break; 549dc008846SHong Zhang if (PetscBTLookup(table0,col)) { 55026fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row+rstart)) nidx_i[isz++] = row+rstart; 551dc008846SHong Zhang break; /* for l = start; l<end ; l++) */ 552632d0f97SHong Zhang } 553dc008846SHong Zhang } 554dc008846SHong Zhang } 555bfc6803cSHong Zhang } 556dc008846SHong Zhang 557dc008846SHong Zhang if (i < is_max - 1) { 558fc70d14eSHong Zhang idx_i += n; /* ptr to input is[i+1] array */ 559bfc6803cSHong Zhang nidx_i += isz; /* ptr to output is[i+1] array */ 560dc008846SHong Zhang } 561fc70d14eSHong Zhang nidx[1+i] = isz; /* size of new is[i] */ 5621ab97a2aSSatish Balay } /* for each is */ 56394bacf5dSBarry Smith ierr = PetscBTDestroy(&table0);CHKERRQ(ierr); 564632d0f97SHong Zhang PetscFunctionReturn(0); 565632d0f97SHong Zhang } 566