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 12632d0f97SHong Zhang #undef __FUNCT__ 13632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ" 1413f74950SBarry Smith PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat C,PetscInt is_max,IS is[],PetscInt ov) 15632d0f97SHong Zhang { 166849ba73SBarry Smith PetscErrorCode ierr; 17b11de2a9SHong Zhang PetscInt i,N=C->cmap->N, bs=C->rmap->bs,M=C->rmap->N,Mbs=M/bs,*nidx,isz,iov; 18b11de2a9SHong Zhang IS *is_new,*is_row; 19b11de2a9SHong Zhang Mat *submats; 20b11de2a9SHong Zhang Mat_MPISBAIJ *c=(Mat_MPISBAIJ*)C->data; 21b11de2a9SHong Zhang Mat_SeqSBAIJ *asub_i; 22b11de2a9SHong Zhang PetscBT table; 23b11de2a9SHong Zhang PetscInt *ai,brow,nz,nis,l,nmax,nstages_local,nstages,max_no,pos; 24b11de2a9SHong Zhang const PetscInt *idx; 254da72fa9SHong Zhang PetscBool flg,*allcolumns,*allrows; 26632d0f97SHong Zhang 27632d0f97SHong Zhang PetscFunctionBegin; 28c910923dSHong Zhang ierr = PetscMalloc(is_max*sizeof(IS),&is_new);CHKERRQ(ierr); 29632d0f97SHong Zhang /* Convert the indices into block format */ 3005d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,is_max,is,is_new);CHKERRQ(ierr); 3126fbe8dcSKarl Rupp if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n"); 32db41cccfSHong Zhang 33b11de2a9SHong Zhang /* ----- previous non-scalable implementation ----- */ 34a0769a71SSatish Balay flg = PETSC_FALSE; 350298fd71SBarry Smith ierr = PetscOptionsHasName(NULL, "-IncreaseOverlap_old", &flg);CHKERRQ(ierr); 367a868f3eSHong Zhang if (flg) { /* previous non-scalable implementation */ 377a868f3eSHong Zhang printf("use previous non-scalable implementation...\n"); 38632d0f97SHong Zhang for (i=0; i<ov; ++i) { 39c910923dSHong Zhang ierr = MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr); 40632d0f97SHong Zhang } 411c4982ebSHong Zhang } else { /* implementation using modified BAIJ routines */ 421c4982ebSHong Zhang 43db41cccfSHong Zhang ierr = PetscMalloc((Mbs+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr); 4453b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&table);CHKERRQ(ierr); /* for column search */ 45*dcca6d9dSJed Brown ierr = PetscMalloc2(is_max+1,&allcolumns,is_max+1,&allrows);CHKERRQ(ierr); 46db41cccfSHong Zhang 47db41cccfSHong Zhang /* Create is_row */ 48db41cccfSHong Zhang ierr = PetscMalloc(is_max*sizeof(IS **),&is_row);CHKERRQ(ierr); 494da72fa9SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,Mbs,0,1,&is_row[0]);CHKERRQ(ierr); 5026fbe8dcSKarl Rupp 511c4982ebSHong Zhang allrows[0] = PETSC_TRUE; 524da72fa9SHong Zhang for (i=1; i<is_max; i++) { 534da72fa9SHong Zhang is_row[i] = is_row[0]; /* reuse is_row[0] */ 544da72fa9SHong Zhang allrows[i] = PETSC_TRUE; 554da72fa9SHong Zhang } 56db41cccfSHong Zhang 57b11de2a9SHong Zhang /* Allocate memory to hold all the submatrices - Modified from MatGetSubMatrices_MPIBAIJ() */ 58db41cccfSHong Zhang ierr = PetscMalloc((is_max+1)*sizeof(Mat),&submats);CHKERRQ(ierr); 59b11de2a9SHong Zhang 60307b7a18SHong Zhang /* Check for special case: each processor gets entire matrix columns */ 61307b7a18SHong Zhang for (i=0; i<is_max; i++) { 62307b7a18SHong Zhang ierr = ISIdentity(is_new[i],&flg);CHKERRQ(ierr); 63307b7a18SHong Zhang ierr = ISGetLocalSize(is_new[i],&isz);CHKERRQ(ierr); 64307b7a18SHong Zhang if (flg && isz == Mbs) { 65307b7a18SHong Zhang allcolumns[i] = PETSC_TRUE; 66307b7a18SHong Zhang } else { 67307b7a18SHong Zhang allcolumns[i] = PETSC_FALSE; 68307b7a18SHong Zhang } 69307b7a18SHong Zhang } 70307b7a18SHong Zhang 71db41cccfSHong Zhang /* Determine the number of stages through which submatrices are done */ 72db41cccfSHong Zhang nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 73db41cccfSHong Zhang if (!nmax) nmax = 1; 74db41cccfSHong Zhang nstages_local = is_max/nmax + ((is_max % nmax) ? 1 : 0); 75db41cccfSHong Zhang 76db41cccfSHong Zhang /* Make sure every processor loops through the nstages */ 77ce94432eSBarry Smith ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 78b11de2a9SHong Zhang 79b11de2a9SHong Zhang for (iov=0; iov<ov; ++iov) { 80b11de2a9SHong Zhang /* 1) Get submats for column search */ 81db41cccfSHong Zhang for (i=0,pos=0; i<nstages; i++) { 82db41cccfSHong Zhang if (pos+nmax <= is_max) max_no = nmax; 83db41cccfSHong Zhang else if (pos == is_max) max_no = 0; 84db41cccfSHong Zhang else max_no = is_max-pos; 857a868f3eSHong Zhang c->ijonly = PETSC_TRUE; 861c4982ebSHong Zhang ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,is_row+pos,is_new+pos,MAT_INITIAL_MATRIX,allrows+pos,allcolumns+pos,submats+pos);CHKERRQ(ierr); 87db41cccfSHong Zhang pos += max_no; 88db41cccfSHong Zhang } 89db41cccfSHong Zhang 90b11de2a9SHong Zhang /* 2) Row search */ 91db41cccfSHong Zhang ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr); 92db41cccfSHong Zhang 93b11de2a9SHong Zhang /* 3) Column search */ 94db41cccfSHong Zhang for (i=0; i<is_max; i++) { 95db41cccfSHong Zhang asub_i = (Mat_SeqSBAIJ*)submats[i]->data; 96db41cccfSHong Zhang ai = asub_i->i;; 97db41cccfSHong Zhang 98db41cccfSHong Zhang /* put is_new obtained from MatIncreaseOverlap_MPIBAIJ() to table */ 99db41cccfSHong Zhang ierr = PetscBTMemzero(Mbs,table);CHKERRQ(ierr); 100db41cccfSHong Zhang 101db41cccfSHong Zhang ierr = ISGetIndices(is_new[i],&idx);CHKERRQ(ierr); 102db41cccfSHong Zhang ierr = ISGetLocalSize(is_new[i],&nis);CHKERRQ(ierr); 103db41cccfSHong Zhang for (l=0; l<nis; l++) { 104db41cccfSHong Zhang ierr = PetscBTSet(table,idx[l]);CHKERRQ(ierr); 105db41cccfSHong Zhang nidx[l] = idx[l]; 106db41cccfSHong Zhang } 107db41cccfSHong Zhang isz = nis; 108db41cccfSHong Zhang 109b11de2a9SHong Zhang /* add column entries to table */ 110db41cccfSHong Zhang for (brow=0; brow<Mbs; brow++) { 111db41cccfSHong Zhang nz = ai[brow+1] - ai[brow]; 112db41cccfSHong Zhang if (nz) { 113db41cccfSHong Zhang if (!PetscBTLookupSet(table,brow)) nidx[isz++] = brow; 114db41cccfSHong Zhang } 115db41cccfSHong Zhang } 116db41cccfSHong Zhang ierr = ISRestoreIndices(is_new[i],&idx);CHKERRQ(ierr); 1176bf464f9SBarry Smith ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr); 118db41cccfSHong Zhang 119db41cccfSHong Zhang /* create updated is_new */ 120db41cccfSHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,is_new+i);CHKERRQ(ierr); 121db41cccfSHong Zhang } 122db41cccfSHong Zhang 123db41cccfSHong Zhang /* Free tmp spaces */ 124db41cccfSHong Zhang for (i=0; i<is_max; i++) { 1256bf464f9SBarry Smith ierr = MatDestroy(&submats[i]);CHKERRQ(ierr); 126db41cccfSHong Zhang } 127db41cccfSHong Zhang } 1284da72fa9SHong Zhang ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr); 12994bacf5dSBarry Smith ierr = PetscBTDestroy(&table);CHKERRQ(ierr); 130db41cccfSHong Zhang ierr = PetscFree(submats);CHKERRQ(ierr); 1316bf464f9SBarry Smith ierr = ISDestroy(&is_row[0]);CHKERRQ(ierr); 132db41cccfSHong Zhang ierr = PetscFree(is_row);CHKERRQ(ierr); 133db41cccfSHong Zhang ierr = PetscFree(nidx);CHKERRQ(ierr); 134b11de2a9SHong Zhang 135db41cccfSHong Zhang } 136db41cccfSHong Zhang 1376bf464f9SBarry Smith for (i=0; i<is_max; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 13805d8c843SHong Zhang ierr = ISExpandIndicesGeneral(N,N,bs,is_max,is_new,is);CHKERRQ(ierr); 139db41cccfSHong Zhang 1406bf464f9SBarry Smith for (i=0; i<is_max; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 141632d0f97SHong Zhang ierr = PetscFree(is_new);CHKERRQ(ierr); 142632d0f97SHong Zhang PetscFunctionReturn(0); 143632d0f97SHong Zhang } 144632d0f97SHong Zhang 1454a69c536SBarry Smith typedef enum {MINE,OTHER} WhoseOwner; 1460472cc68SHong Zhang /* data1, odata1 and odata2 are packed in the format (for communication): 147a2a9f22aSHong Zhang data[0] = is_max, no of is 148a2a9f22aSHong Zhang data[1] = size of is[0] 149a2a9f22aSHong Zhang ... 150a2a9f22aSHong Zhang data[is_max] = size of is[is_max-1] 151a2a9f22aSHong Zhang data[is_max + 1] = data(is[0]) 152a2a9f22aSHong Zhang ... 153a2a9f22aSHong Zhang data[is_max+1+sum(size of is[k]), k=0,...,i-1] = data(is[i]) 154a2a9f22aSHong Zhang ... 1550472cc68SHong Zhang data2 is packed in the format (for creating output is[]): 1560472cc68SHong Zhang data[0] = is_max, no of is 1570472cc68SHong Zhang data[1] = size of is[0] 1580472cc68SHong Zhang ... 1590472cc68SHong Zhang data[is_max] = size of is[is_max-1] 1600472cc68SHong Zhang data[is_max + 1] = data(is[0]) 1610472cc68SHong Zhang ... 1620472cc68SHong Zhang data[is_max + 1 + Mbs*i) = data(is[i]) 1630472cc68SHong Zhang ... 164a2a9f22aSHong Zhang */ 165632d0f97SHong Zhang #undef __FUNCT__ 166632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Once" 16713f74950SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat C,PetscInt is_max,IS is[]) 168632d0f97SHong Zhang { 169632d0f97SHong Zhang Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 1706849ba73SBarry Smith PetscErrorCode ierr; 171adcec1e5SJed Brown PetscMPIInt size,rank,tag1,tag2,*len_s,nrqr,nrqs,*id_r1,*len_r1,flag,len,*iwork; 1725d0c19d7SBarry Smith const PetscInt *idx_i; 17326fbe8dcSKarl Rupp PetscInt idx,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i; 17426fbe8dcSKarl Rupp PetscInt Mbs,i,j,k,*odata1,*odata2; 17526fbe8dcSKarl Rupp PetscInt proc_id,**odata2_ptr,*ctable=0,*btable,len_max,len_est; 176adcec1e5SJed Brown PetscInt proc_end=0,len_unused,nodata2; 17713f74950SBarry Smith PetscInt ois_max; /* max no of is[] in each of processor */ 178bfc6803cSHong Zhang char *t_p; 179632d0f97SHong Zhang MPI_Comm comm; 180e8527bf2SHong Zhang MPI_Request *s_waits1,*s_waits2,r_req; 181632d0f97SHong Zhang MPI_Status *s_status,r_status; 18245f43ab7SHong Zhang PetscBT *table; /* mark indices of this processor's is[] */ 183fc70d14eSHong Zhang PetscBT table_i; 184bfc6803cSHong Zhang PetscBT otable; /* mark indices of other processors' is[] */ 185d0f46423SBarry Smith PetscInt bs=C->rmap->bs,Bn = c->B->cmap->n,Bnbs = Bn/bs,*Bowners; 18610eea884SHong Zhang IS garray_local,garray_gl; 1875483b11dSHong Zhang 188632d0f97SHong Zhang PetscFunctionBegin; 189ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 190632d0f97SHong Zhang size = c->size; 191632d0f97SHong Zhang rank = c->rank; 192632d0f97SHong Zhang Mbs = c->Mbs; 193632d0f97SHong Zhang 194c910923dSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 195c910923dSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 196c910923dSHong Zhang 197430a0127SHong Zhang /* create tables used in 198430a0127SHong Zhang step 1: table[i] - mark c->garray of proc [i] 19945f43ab7SHong Zhang step 3: table[i] - mark indices of is[i] when whose=MINE 200430a0127SHong Zhang table[0] - mark incideces of is[] when whose=OTHER */ 201430a0127SHong Zhang len = PetscMax(is_max, size);CHKERRQ(ierr); 202*dcca6d9dSJed Brown ierr = PetscMalloc2(len,&table,(Mbs/PETSC_BITS_PER_BYTE+1)*len,&t_p);CHKERRQ(ierr); 203430a0127SHong Zhang for (i=0; i<len; i++) { 204430a0127SHong Zhang table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 205430a0127SHong Zhang } 206430a0127SHong Zhang 20713f74950SBarry Smith ierr = MPI_Allreduce(&is_max,&ois_max,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 20810eea884SHong Zhang 2095483b11dSHong Zhang /* 1. Send this processor's is[] to other processors */ 2105483b11dSHong Zhang /*---------------------------------------------------*/ 211e8527bf2SHong Zhang /* allocate spaces */ 21213f74950SBarry Smith ierr = PetscMalloc(is_max*sizeof(PetscInt),&n);CHKERRQ(ierr); 2135483b11dSHong Zhang len = 0; 214c910923dSHong Zhang for (i=0; i<is_max; i++) { 215632d0f97SHong Zhang ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 216632d0f97SHong Zhang len += n[i]; 217632d0f97SHong Zhang } 218958c9bccSBarry Smith if (!len) { 2195483b11dSHong Zhang is_max = 0; 2205483b11dSHong Zhang } else { 2215483b11dSHong Zhang len += 1 + is_max; /* max length of data1 for one processor */ 2225483b11dSHong Zhang } 2235483b11dSHong Zhang 224430a0127SHong Zhang 22513f74950SBarry Smith ierr = PetscMalloc((size*len+1)*sizeof(PetscInt),&data1);CHKERRQ(ierr); 22613f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt*),&data1_start);CHKERRQ(ierr); 2275483b11dSHong Zhang for (i=0; i<size; i++) data1_start[i] = data1 + i*len; 2285483b11dSHong Zhang 229*dcca6d9dSJed Brown ierr = PetscMalloc4(size,&len_s,size,&btable,size,&iwork,size+1,&Bowners);CHKERRQ(ierr); 230e8527bf2SHong Zhang 23176f244e2SHong Zhang /* gather c->garray from all processors */ 23270b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,Bnbs,c->garray,PETSC_COPY_VALUES,&garray_local);CHKERRQ(ierr); 23376f244e2SHong Zhang ierr = ISAllGather(garray_local, &garray_gl);CHKERRQ(ierr); 2346bf464f9SBarry Smith ierr = ISDestroy(&garray_local);CHKERRQ(ierr); 235a7cc72afSBarry Smith ierr = MPI_Allgather(&Bnbs,1,MPIU_INT,Bowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 23626fbe8dcSKarl Rupp 23776f244e2SHong Zhang Bowners[0] = 0; 23876f244e2SHong Zhang for (i=0; i<size; i++) Bowners[i+1] += Bowners[i]; 23976f244e2SHong Zhang 2405483b11dSHong Zhang if (is_max) { 241430a0127SHong Zhang /* hash table ctable which maps c->row to proc_id) */ 24213f74950SBarry Smith ierr = PetscMalloc(Mbs*sizeof(PetscInt),&ctable);CHKERRQ(ierr); 2435483b11dSHong Zhang for (proc_id=0,j=0; proc_id<size; proc_id++) { 24426fbe8dcSKarl Rupp for (; j<C->rmap->range[proc_id+1]/bs; j++) ctable[j] = proc_id; 2455483b11dSHong Zhang } 2465483b11dSHong Zhang 247430a0127SHong Zhang /* hash tables marking c->garray */ 24822d28d08SBarry Smith ierr = ISGetIndices(garray_gl,&idx_i);CHKERRQ(ierr); 2495483b11dSHong Zhang for (i=0; i<size; i++) { 250430a0127SHong Zhang table_i = table[i]; 251430a0127SHong Zhang ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); 252430a0127SHong Zhang for (j = Bowners[i]; j<Bowners[i+1]; j++) { /* go through B cols of proc[i]*/ 253430a0127SHong Zhang ierr = PetscBTSet(table_i,idx_i[j]);CHKERRQ(ierr); 2545483b11dSHong Zhang } 2555483b11dSHong Zhang } 25610eea884SHong Zhang ierr = ISRestoreIndices(garray_gl,&idx_i);CHKERRQ(ierr); 2575483b11dSHong Zhang } /* if (is_max) */ 2586bf464f9SBarry Smith ierr = ISDestroy(&garray_gl);CHKERRQ(ierr); 2595483b11dSHong Zhang 2605483b11dSHong Zhang /* evaluate communication - mesg to who, length, and buffer space */ 261e8527bf2SHong Zhang for (i=0; i<size; i++) len_s[i] = 0; 2625483b11dSHong Zhang 2635483b11dSHong Zhang /* header of data1 */ 2645483b11dSHong Zhang for (proc_id=0; proc_id<size; proc_id++) { 2655483b11dSHong Zhang iwork[proc_id] = 0; 2665483b11dSHong Zhang *data1_start[proc_id] = is_max; 2675483b11dSHong Zhang data1_start[proc_id]++; 2685483b11dSHong Zhang for (j=0; j<is_max; j++) { 2695483b11dSHong Zhang if (proc_id == rank) { 2705483b11dSHong Zhang *data1_start[proc_id] = n[j]; 2715483b11dSHong Zhang } else { 2725483b11dSHong Zhang *data1_start[proc_id] = 0; 2735483b11dSHong Zhang } 2745483b11dSHong Zhang data1_start[proc_id]++; 2755483b11dSHong Zhang } 2765483b11dSHong Zhang } 2775483b11dSHong Zhang 2785483b11dSHong Zhang for (i=0; i<is_max; i++) { 2795483b11dSHong Zhang ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr); 2805483b11dSHong Zhang for (j=0; j<n[i]; j++) { 2815483b11dSHong Zhang idx = idx_i[j]; 2825483b11dSHong Zhang *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */ 2835483b11dSHong Zhang proc_end = ctable[idx]; 2845483b11dSHong Zhang for (proc_id=0; proc_id<=proc_end; proc_id++) { /* for others to process */ 285e8527bf2SHong Zhang if (proc_id == rank) continue; /* done before this loop */ 28626fbe8dcSKarl Rupp if (proc_id < proc_end && !PetscBTLookup(table[proc_id],idx)) continue; /* no need for sending idx to [proc_id] */ 2875483b11dSHong Zhang *data1_start[proc_id] = idx; data1_start[proc_id]++; 2885483b11dSHong Zhang len_s[proc_id]++; 2895483b11dSHong Zhang } 2905483b11dSHong Zhang } 2915483b11dSHong Zhang /* update header data */ 2922cfbe0a4SHong Zhang for (proc_id=0; proc_id<size; proc_id++) { 2935483b11dSHong Zhang if (proc_id== rank) continue; 2945483b11dSHong Zhang *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id]; 2955483b11dSHong Zhang iwork[proc_id] = len_s[proc_id]; 2965483b11dSHong Zhang } 2975483b11dSHong Zhang ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr); 298e8527bf2SHong Zhang } 2995483b11dSHong Zhang 3005483b11dSHong Zhang nrqs = 0; nrqr = 0; 3015483b11dSHong Zhang for (i=0; i<size; i++) { 3025483b11dSHong Zhang data1_start[i] = data1 + i*len; 3035483b11dSHong Zhang if (len_s[i]) { 3045483b11dSHong Zhang nrqs++; 3055483b11dSHong Zhang len_s[i] += 1 + is_max; /* add no. of header msg */ 3065483b11dSHong Zhang } 3075483b11dSHong Zhang } 3085483b11dSHong Zhang 3095483b11dSHong Zhang for (i=0; i<is_max; i++) { 3106bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 311632d0f97SHong Zhang } 312bfc6803cSHong Zhang ierr = PetscFree(n);CHKERRQ(ierr); 31305b42c5fSBarry Smith ierr = PetscFree(ctable);CHKERRQ(ierr); 3145483b11dSHong Zhang 3155483b11dSHong Zhang /* Determine the number of messages to expect, their lengths, from from-ids */ 3160298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrqr);CHKERRQ(ierr); 3175483b11dSHong Zhang ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);CHKERRQ(ierr); 318632d0f97SHong Zhang 319632d0f97SHong Zhang /* Now post the sends */ 320*dcca6d9dSJed Brown ierr = PetscMalloc2(size,&s_waits1,size,&s_waits2);CHKERRQ(ierr); 321632d0f97SHong Zhang k = 0; 3225483b11dSHong Zhang for (proc_id=0; proc_id<size; proc_id++) { /* send data1 to processor [proc_id] */ 3235483b11dSHong Zhang if (len_s[proc_id]) { 32413f74950SBarry Smith ierr = MPI_Isend(data1_start[proc_id],len_s[proc_id],MPIU_INT,proc_id,tag1,comm,s_waits1+k);CHKERRQ(ierr); 325632d0f97SHong Zhang k++; 326632d0f97SHong Zhang } 327632d0f97SHong Zhang } 328632d0f97SHong Zhang 32945f43ab7SHong Zhang /* 2. Receive other's is[] and process. Then send back */ 330bfc6803cSHong Zhang /*-----------------------------------------------------*/ 3315483b11dSHong Zhang len = 0; 3325483b11dSHong Zhang for (i=0; i<nrqr; i++) { 3335483b11dSHong Zhang if (len_r1[i] > len) len = len_r1[i]; 3345483b11dSHong Zhang } 33545f43ab7SHong Zhang ierr = PetscFree(len_r1);CHKERRQ(ierr); 33645f43ab7SHong Zhang ierr = PetscFree(id_r1);CHKERRQ(ierr); 33745f43ab7SHong Zhang 33826fbe8dcSKarl Rupp for (proc_id=0; proc_id<size; proc_id++) len_s[proc_id] = iwork[proc_id] = 0; 33945f43ab7SHong Zhang 34013f74950SBarry Smith ierr = PetscMalloc((len+1)*sizeof(PetscInt),&odata1);CHKERRQ(ierr); 34113f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt**),&odata2_ptr);CHKERRQ(ierr); 34253b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&otable);CHKERRQ(ierr); 34310eea884SHong Zhang 34410eea884SHong Zhang len_max = ois_max*(Mbs+1); /* max space storing all is[] for each receive */ 345240e5896SHong Zhang len_est = 2*len_max; /* estimated space of storing is[] for all receiving messages */ 34613f74950SBarry Smith ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr); 34710eea884SHong Zhang nodata2 = 0; /* nodata2+1: num of PetscMalloc(,&odata2_ptr[]) called */ 34826fbe8dcSKarl Rupp 349240e5896SHong Zhang odata2_ptr[nodata2] = odata2; 35026fbe8dcSKarl Rupp 35110eea884SHong Zhang len_unused = len_est; /* unused space in the array odata2_ptr[nodata2]-- needs to be >= len_max */ 35210eea884SHong Zhang 353632d0f97SHong Zhang k = 0; 3545483b11dSHong Zhang while (k < nrqr) { 355632d0f97SHong Zhang /* Receive messages */ 356bfc6803cSHong Zhang ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);CHKERRQ(ierr); 357632d0f97SHong Zhang if (flag) { 35813f74950SBarry Smith ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRQ(ierr); 359632d0f97SHong Zhang proc_id = r_status.MPI_SOURCE; 36013f74950SBarry Smith ierr = MPI_Irecv(odata1,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr); 3615483b11dSHong Zhang ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr); 362632d0f97SHong Zhang 363fc70d14eSHong Zhang /* Process messages */ 364240e5896SHong Zhang /* make sure there is enough unused space in odata2 array */ 36510eea884SHong Zhang if (len_unused < len_max) { /* allocate more space for odata2 */ 36613f74950SBarry Smith ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr); 36726fbe8dcSKarl Rupp 368240e5896SHong Zhang odata2_ptr[++nodata2] = odata2; 36926fbe8dcSKarl Rupp 37010eea884SHong Zhang len_unused = len_est; 37110eea884SHong Zhang } 37210eea884SHong Zhang 37310eea884SHong Zhang ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,odata2,&otable);CHKERRQ(ierr); 374a2a9f22aSHong Zhang len = 1 + odata2[0]; 37526fbe8dcSKarl Rupp for (i=0; i<odata2[0]; i++) len += odata2[1 + i]; 376632d0f97SHong Zhang 377632d0f97SHong Zhang /* Send messages back */ 37813f74950SBarry Smith ierr = MPI_Isend(odata2,len,MPIU_INT,proc_id,tag2,comm,s_waits2+k);CHKERRQ(ierr); 379fc70d14eSHong Zhang k++; 38010eea884SHong Zhang odata2 += len; 38110eea884SHong Zhang len_unused -= len; 38245f43ab7SHong Zhang len_s[proc_id] = len; /* num of messages sending back to [proc_id] by this proc */ 383632d0f97SHong Zhang } 3845483b11dSHong Zhang } 3855483b11dSHong Zhang ierr = PetscFree(odata1);CHKERRQ(ierr); 38694bacf5dSBarry Smith ierr = PetscBTDestroy(&otable);CHKERRQ(ierr); 387632d0f97SHong Zhang 38845f43ab7SHong Zhang /* 3. Do local work on this processor's is[] */ 38945f43ab7SHong Zhang /*-------------------------------------------*/ 39045f43ab7SHong Zhang /* make sure there is enough unused space in odata2(=data) array */ 39145f43ab7SHong Zhang len_max = is_max*(Mbs+1); /* max space storing all is[] for this processor */ 39210eea884SHong Zhang if (len_unused < len_max) { /* allocate more space for odata2 */ 39313f74950SBarry Smith ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr); 39426fbe8dcSKarl Rupp 395240e5896SHong Zhang odata2_ptr[++nodata2] = odata2; 39610eea884SHong Zhang } 397bfc6803cSHong Zhang 398240e5896SHong Zhang data = odata2; 39945f43ab7SHong Zhang ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,data,table);CHKERRQ(ierr); 40045f43ab7SHong Zhang ierr = PetscFree(data1_start);CHKERRQ(ierr); 40145f43ab7SHong Zhang 40245f43ab7SHong Zhang /* 4. Receive work done on other processors, then merge */ 40345f43ab7SHong Zhang /*------------------------------------------------------*/ 40445f43ab7SHong Zhang /* get max number of messages that this processor expects to recv */ 405adcec1e5SJed Brown ierr = MPI_Allreduce(len_s,iwork,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 40613f74950SBarry Smith ierr = PetscMalloc((iwork[rank]+1)*sizeof(PetscInt),&data2);CHKERRQ(ierr); 40774ed9c26SBarry Smith ierr = PetscFree4(len_s,btable,iwork,Bowners);CHKERRQ(ierr); 40845f43ab7SHong Zhang 409632d0f97SHong Zhang k = 0; 4105483b11dSHong Zhang while (k < nrqs) { 411632d0f97SHong Zhang /* Receive messages */ 41222d28d08SBarry Smith ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status);CHKERRQ(ierr); 413632d0f97SHong Zhang if (flag) { 41413f74950SBarry Smith ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRQ(ierr); 41526fbe8dcSKarl Rupp 416632d0f97SHong Zhang proc_id = r_status.MPI_SOURCE; 41726fbe8dcSKarl Rupp 41813f74950SBarry Smith ierr = MPI_Irecv(data2,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr); 4195483b11dSHong Zhang ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr); 4205483b11dSHong Zhang if (len > 1+is_max) { /* Add data2 into data */ 4210472cc68SHong Zhang data2_i = data2 + 1 + is_max; 422fc70d14eSHong Zhang for (i=0; i<is_max; i++) { 423fc70d14eSHong Zhang table_i = table[i]; 424bfc6803cSHong Zhang data_i = data + 1 + is_max + Mbs*i; 425bfc6803cSHong Zhang isz = data[1+i]; 4260472cc68SHong Zhang for (j=0; j<data2[1+i]; j++) { 4270472cc68SHong Zhang col = data2_i[j]; 42826fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,col)) data_i[isz++] = col; 429fc70d14eSHong Zhang } 430bfc6803cSHong Zhang data[1+i] = isz; 4310472cc68SHong Zhang if (i < is_max - 1) data2_i += data2[1+i]; 432fc70d14eSHong Zhang } 4335483b11dSHong Zhang } 434632d0f97SHong Zhang k++; 435632d0f97SHong Zhang } 4365483b11dSHong Zhang } 43745f43ab7SHong Zhang ierr = PetscFree(data2);CHKERRQ(ierr); 43874ed9c26SBarry Smith ierr = PetscFree2(table,t_p);CHKERRQ(ierr); 4395483b11dSHong Zhang 4405483b11dSHong Zhang /* phase 1 sends are complete */ 4415483b11dSHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 4420c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 4435483b11dSHong Zhang ierr = PetscFree(data1);CHKERRQ(ierr); 4445483b11dSHong Zhang 445240e5896SHong Zhang /* phase 2 sends are complete */ 4460c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status);CHKERRQ(ierr);} 44774ed9c26SBarry Smith ierr = PetscFree2(s_waits1,s_waits2);CHKERRQ(ierr); 44845f43ab7SHong Zhang ierr = PetscFree(s_status);CHKERRQ(ierr); 449632d0f97SHong Zhang 450c910923dSHong Zhang /* 5. Create new is[] */ 451c910923dSHong Zhang /*--------------------*/ 452c910923dSHong Zhang for (i=0; i<is_max; i++) { 453bfc6803cSHong Zhang data_i = data + 1 + is_max + Mbs*i; 45470b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 455632d0f97SHong Zhang } 45645f43ab7SHong Zhang for (k=0; k<=nodata2; k++) { 45745f43ab7SHong Zhang ierr = PetscFree(odata2_ptr[k]);CHKERRQ(ierr); 45845f43ab7SHong Zhang } 45945f43ab7SHong Zhang ierr = PetscFree(odata2_ptr);CHKERRQ(ierr); 460632d0f97SHong Zhang PetscFunctionReturn(0); 461632d0f97SHong Zhang } 462632d0f97SHong Zhang 463632d0f97SHong Zhang #undef __FUNCT__ 464632d0f97SHong Zhang #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Local" 465632d0f97SHong Zhang /* 466dc008846SHong Zhang MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do 467632d0f97SHong Zhang the work on the local processor. 468632d0f97SHong Zhang 469632d0f97SHong Zhang Inputs: 470632d0f97SHong Zhang C - MAT_MPISBAIJ; 471bfc6803cSHong Zhang data - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format. 472bfc6803cSHong Zhang whose - whose is[] to be processed, 473bfc6803cSHong Zhang MINE: this processor's is[] 474bfc6803cSHong Zhang OTHER: other processor's is[] 475632d0f97SHong Zhang Output: 47610eea884SHong Zhang nidx - whose = MINE: 4770472cc68SHong Zhang holds input and newly found indices in the same format as data 4780472cc68SHong Zhang whose = OTHER: 4790472cc68SHong Zhang only holds the newly found indices 4800472cc68SHong Zhang table - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'. 481632d0f97SHong Zhang */ 48276f244e2SHong Zhang /* Would computation be reduced by swapping the loop 'for each is' and 'for each row'? */ 48313f74950SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat C,PetscInt *data,PetscInt whose,PetscInt *nidx,PetscBT *table) 484632d0f97SHong Zhang { 485632d0f97SHong Zhang Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 486dc008846SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(c->A)->data; 487dc008846SHong Zhang Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(c->B)->data; 488dfbe8321SBarry Smith PetscErrorCode ierr; 48913f74950SBarry Smith PetscInt row,mbs,Mbs,*nidx_i,col,col_max,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l; 49013f74950SBarry Smith PetscInt a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n; 491bfc6803cSHong Zhang PetscBT table0; /* mark the indices of input is[] for look up */ 492bfc6803cSHong Zhang PetscBT table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */ 493632d0f97SHong Zhang 494632d0f97SHong Zhang PetscFunctionBegin; 49531f99336SHong Zhang Mbs = c->Mbs; mbs = a->mbs; 496dc008846SHong Zhang ai = a->i; aj = a->j; 497dc008846SHong Zhang bi = b->i; bj = b->j; 498632d0f97SHong Zhang garray = c->garray; 499899cda47SBarry Smith rstart = c->rstartbs; 500dc008846SHong Zhang is_max = data[0]; 501632d0f97SHong Zhang 50253b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&table0);CHKERRQ(ierr); 503fc70d14eSHong Zhang 504fc70d14eSHong Zhang nidx[0] = is_max; 505fc70d14eSHong Zhang idx_i = data + is_max + 1; /* ptr to input is[0] array */ 506bfc6803cSHong Zhang nidx_i = nidx + is_max + 1; /* ptr to output is[0] array */ 507dc008846SHong Zhang for (i=0; i<is_max; i++) { /* for each is */ 508dc008846SHong Zhang isz = 0; 509fc70d14eSHong Zhang n = data[1+i]; /* size of input is[i] */ 510dc008846SHong Zhang 51176f244e2SHong Zhang /* initialize and set table_i(mark idx and nidx) and table0(only mark idx) */ 512bfc6803cSHong Zhang if (whose == MINE) { /* process this processor's is[] */ 513bfc6803cSHong Zhang table_i = table[i]; 5140472cc68SHong Zhang nidx_i = nidx + 1+ is_max + Mbs*i; 515bfc6803cSHong Zhang } else { /* process other processor's is[] - only use one temp table */ 516430a0127SHong Zhang table_i = table[0]; 517bfc6803cSHong Zhang } 518bfc6803cSHong Zhang ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); 519bfc6803cSHong Zhang ierr = PetscBTMemzero(Mbs,table0);CHKERRQ(ierr); 52076f244e2SHong Zhang if (n==0) { 52176f244e2SHong Zhang nidx[1+i] = 0; /* size of new is[i] */ 52276f244e2SHong Zhang continue; 52376f244e2SHong Zhang } 52476f244e2SHong Zhang 52576f244e2SHong Zhang isz0 = 0; col_max = 0; 526dc008846SHong Zhang for (j=0; j<n; j++) { 527dc008846SHong Zhang col = idx_i[j]; 528e32f2f54SBarry Smith if (col >= Mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index col %D >= Mbs %D",col,Mbs); 529bfc6803cSHong Zhang if (!PetscBTLookupSet(table_i,col)) { 530bfc6803cSHong Zhang ierr = PetscBTSet(table0,col);CHKERRQ(ierr); 53126fbe8dcSKarl Rupp if (whose == MINE) nidx_i[isz0] = col; 532dbe03f88SHong Zhang if (col_max < col) col_max = col; 533bfc6803cSHong Zhang isz0++; 534bfc6803cSHong Zhang } 535632d0f97SHong Zhang } 536dc008846SHong Zhang 53726fbe8dcSKarl Rupp if (whose == MINE) isz = isz0; 538fc70d14eSHong Zhang k = 0; /* no. of indices from input is[i] that have been examined */ 539dc008846SHong Zhang for (row=0; row<mbs; row++) { 540dc008846SHong Zhang a_start = ai[row]; a_end = ai[row+1]; 541dc008846SHong Zhang b_start = bi[row]; b_end = bi[row+1]; 5420472cc68SHong Zhang if (PetscBTLookup(table0,row+rstart)) { /* row is on input is[i]: 5430472cc68SHong Zhang do row search: collect all col in this row */ 544dc008846SHong Zhang for (l = a_start; l<a_end ; l++) { /* Amat */ 545dc008846SHong Zhang col = aj[l] + rstart; 54626fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,col)) nidx_i[isz++] = col; 547dc008846SHong Zhang } 548dc008846SHong Zhang for (l = b_start; l<b_end ; l++) { /* Bmat */ 549dc008846SHong Zhang col = garray[bj[l]]; 55026fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,col)) nidx_i[isz++] = col; 551dc008846SHong Zhang } 552dc008846SHong Zhang k++; 553dc008846SHong Zhang if (k >= isz0) break; /* for (row=0; row<mbs; row++) */ 5540472cc68SHong Zhang } else { /* row is not on input is[i]: 5550472cc68SHong Zhang do col serach: add row onto nidx_i if there is a col in nidx_i */ 556dc008846SHong Zhang for (l = a_start; l<a_end; l++) { /* Amat */ 557dc008846SHong Zhang col = aj[l] + rstart; 55876f244e2SHong Zhang if (col > col_max) break; 559dc008846SHong Zhang if (PetscBTLookup(table0,col)) { 56026fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row+rstart)) nidx_i[isz++] = row+rstart; 561dc008846SHong Zhang break; /* for l = start; l<end ; l++) */ 562632d0f97SHong Zhang } 563632d0f97SHong Zhang } 564dc008846SHong Zhang for (l = b_start; l<b_end; l++) { /* Bmat */ 565dc008846SHong Zhang col = garray[bj[l]]; 56676f244e2SHong Zhang if (col > col_max) break; 567dc008846SHong Zhang if (PetscBTLookup(table0,col)) { 56826fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row+rstart)) nidx_i[isz++] = row+rstart; 569dc008846SHong Zhang break; /* for l = start; l<end ; l++) */ 570632d0f97SHong Zhang } 571dc008846SHong Zhang } 572dc008846SHong Zhang } 573bfc6803cSHong Zhang } 574dc008846SHong Zhang 575dc008846SHong Zhang if (i < is_max - 1) { 576fc70d14eSHong Zhang idx_i += n; /* ptr to input is[i+1] array */ 577bfc6803cSHong Zhang nidx_i += isz; /* ptr to output is[i+1] array */ 578dc008846SHong Zhang } 579fc70d14eSHong Zhang nidx[1+i] = isz; /* size of new is[i] */ 5801ab97a2aSSatish Balay } /* for each is */ 58194bacf5dSBarry Smith ierr = PetscBTDestroy(&table0);CHKERRQ(ierr); 582632d0f97SHong Zhang PetscFunctionReturn(0); 583632d0f97SHong Zhang } 584632d0f97SHong Zhang 585632d0f97SHong Zhang 586