123bdbc58SBarry Smith 2d5b485abSSatish Balay /* 3d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 4d5b485abSSatish Balay and to find submatrices that were shared across processors. 5d5b485abSSatish Balay */ 6c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h> 7c6db04a5SJed Brown #include <petscbt.h> 8d5b485abSSatish Balay 9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char **,PetscInt*,PetscInt**); 10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt **,PetscInt**,PetscInt*); 1109573ac7SBarry Smith extern PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 1209573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 13d5b485abSSatish Balay 144a2ae208SSatish Balay #undef __FUNCT__ 154a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ" 16b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov) 17d5b485abSSatish Balay { 186849ba73SBarry Smith PetscErrorCode ierr; 19d0f46423SBarry Smith PetscInt i,N=C->cmap->N, bs=C->rmap->bs; 2036f4e84dSSatish Balay IS *is_new; 2136f4e84dSSatish Balay 223a40ed3dSBarry Smith PetscFunctionBegin; 2382502324SSatish Balay ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr); 24df36731bSBarry Smith /* Convert the indices into block format */ 2505d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr); 26e32f2f54SBarry Smith if (ov < 0){ SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");} 27d5b485abSSatish Balay for (i=0; i<ov; ++i) { 2836f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 29d5b485abSSatish Balay } 306bf464f9SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 3105d8c843SHong Zhang ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr); 326bf464f9SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 33606d414cSSatish Balay ierr = PetscFree(is_new);CHKERRQ(ierr); 343a40ed3dSBarry Smith PetscFunctionReturn(0); 35d5b485abSSatish Balay } 36d5b485abSSatish Balay 37d5b485abSSatish Balay /* 38d5b485abSSatish Balay Sample message format: 39d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 400e9b0e7eSHong Zhang to index sets is[1], is[5] 41d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 42d5b485abSSatish Balay ----------- 43d5b485abSSatish Balay mesg [1] = 1 => is[1] 44d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 45d5b485abSSatish Balay ----------- 46d5b485abSSatish Balay mesg [5] = 5 => is[5] 47d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 48d5b485abSSatish Balay ----------- 49d5b485abSSatish Balay mesg [7] 500e9b0e7eSHong Zhang mesg [n] data(is[1]) 51d5b485abSSatish Balay ----------- 52d5b485abSSatish Balay mesg[n+1] 53d5b485abSSatish Balay mesg[m] data(is[5]) 54d5b485abSSatish Balay ----------- 55d5b485abSSatish Balay 56d5b485abSSatish Balay Notes: 57d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 58d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 59d5b485abSSatish Balay */ 604a2ae208SSatish Balay #undef __FUNCT__ 614a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once" 62db41cccfSHong Zhang PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[]) 63d5b485abSSatish Balay { 64df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 655d0c19d7SBarry Smith const PetscInt **idx,*idx_i; 6624c049a4SHong Zhang PetscInt *n,*w3,*w4,**data,len; 676849ba73SBarry Smith PetscErrorCode ierr; 68b24ad042SBarry Smith PetscMPIInt size,rank,tag1,tag2,*w2,*w1,nrqr; 69245d216aSHong Zhang PetscInt Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr; 708f9f447aSHong Zhang PetscInt *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p; 71b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1,*onodes2,*olengths2; 72f1af5d2fSBarry Smith PetscBT *table; 73d5b485abSSatish Balay MPI_Comm comm; 74d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 75d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 768f9f447aSHong Zhang char *t_p; 77d5b485abSSatish Balay 783a40ed3dSBarry Smith PetscFunctionBegin; 797adad957SLisandro Dalcin comm = ((PetscObject)C)->comm; 80d5b485abSSatish Balay size = c->size; 81d5b485abSSatish Balay rank = c->rank; 82df36731bSBarry Smith Mbs = c->Mbs; 83d5b485abSSatish Balay 84c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 85c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 86c7dd2924SSatish Balay 8724c049a4SHong Zhang ierr = PetscMalloc2(imax+1,const PetscInt*,&idx,imax,PetscInt,&n);CHKERRQ(ierr); 8824c049a4SHong Zhang 89d5b485abSSatish Balay for (i=0; i<imax; i++) { 90d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 91b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 92d5b485abSSatish Balay } 93d5b485abSSatish Balay 94d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 95d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 9623bdbc58SBarry Smith ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr); 9723bdbc58SBarry Smith ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 9823bdbc58SBarry Smith ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 9923bdbc58SBarry Smith ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr); 100d5b485abSSatish Balay for (i=0; i<imax; i++) { 101b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 102d5b485abSSatish Balay idx_i = idx[i]; 103d5b485abSSatish Balay len = n[i]; 104d5b485abSSatish Balay for (j=0; j<len; j++) { 105d5b485abSSatish Balay row = idx_i[j]; 1066b41c64dSBarry Smith if (row < 0) { 107e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 1086b41c64dSBarry Smith } 109ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 110d5b485abSSatish Balay w4[proc]++; 111d5b485abSSatish Balay } 112d5b485abSSatish Balay for (j=0; j<size; j++){ 113d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 114d5b485abSSatish Balay } 115d5b485abSSatish Balay } 116d5b485abSSatish Balay 117d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 118d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1190e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 120d5b485abSSatish Balay w3[rank] = 0; 121d5b485abSSatish Balay for (i=0; i<size; i++) { 122d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 123d5b485abSSatish Balay } 124d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 125b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); 126d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 127d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 128d5b485abSSatish Balay } 129d5b485abSSatish Balay 130d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 131d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 132d5b485abSSatish Balay j = pa[i]; 133d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 134d5b485abSSatish Balay msz += w1[j]; 135d5b485abSSatish Balay } 136d5b485abSSatish Balay 137c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 138a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 139c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 140d5b485abSSatish Balay 141c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 142c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 143d5b485abSSatish Balay 144d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 14523bdbc58SBarry Smith ierr = PetscMalloc4(size,PetscInt*,&outdat,size,PetscInt*,&ptr,msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr); 14623bdbc58SBarry Smith ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr); 14723bdbc58SBarry Smith ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 148d5b485abSSatish Balay { 149b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 150d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 151d5b485abSSatish Balay j = pa[i]; 152d5b485abSSatish Balay iptr += ict; 153d5b485abSSatish Balay outdat[j] = iptr; 154d5b485abSSatish Balay ict = w1[j]; 155d5b485abSSatish Balay } 156d5b485abSSatish Balay } 157d5b485abSSatish Balay 158d5b485abSSatish Balay /* Form the outgoing messages */ 159d5b485abSSatish Balay /*plug in the headers*/ 160d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 161d5b485abSSatish Balay j = pa[i]; 162d5b485abSSatish Balay outdat[j][0] = 0; 163b24ad042SBarry Smith ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 164d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 165d5b485abSSatish Balay } 166d5b485abSSatish Balay 167d5b485abSSatish Balay /* Memory for doing local proc's work*/ 168d5b485abSSatish Balay { 1698f9f447aSHong Zhang ierr = PetscMalloc5(imax,PetscBT,&table, imax,PetscInt*,&data, imax,PetscInt,&isz, 1708f9f447aSHong Zhang Mbs*imax,PetscInt,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,char,&t_p);CHKERRQ(ierr); 1718f9f447aSHong Zhang ierr = PetscMemzero(table,imax*sizeof(PetscBT));CHKERRQ(ierr); 1728f9f447aSHong Zhang ierr = PetscMemzero(data,imax*sizeof(PetscInt*));CHKERRQ(ierr); 1738f9f447aSHong Zhang ierr = PetscMemzero(isz,imax*sizeof(PetscInt));CHKERRQ(ierr); 1748f9f447aSHong Zhang ierr = PetscMemzero(d_p,Mbs*imax*sizeof(PetscInt));CHKERRQ(ierr); 1758f9f447aSHong Zhang ierr = PetscMemzero(t_p,(Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char));CHKERRQ(ierr); 176d5b485abSSatish Balay 177bbd702dbSSatish Balay for (i=0; i<imax; i++) { 178b6410449SSatish Balay table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 179bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 180d5b485abSSatish Balay } 181d5b485abSSatish Balay } 182d5b485abSSatish Balay 183d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 184d5b485abSSatish Balay { 185b24ad042SBarry Smith PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 186f1af5d2fSBarry Smith PetscBT table_i; 187d5b485abSSatish Balay 188d5b485abSSatish Balay for (i=0; i<imax; i++) { 189b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 190d5b485abSSatish Balay n_i = n[i]; 191d5b485abSSatish Balay table_i = table[i]; 192d5b485abSSatish Balay idx_i = idx[i]; 193d5b485abSSatish Balay data_i = data[i]; 194d5b485abSSatish Balay isz_i = isz[i]; 195d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 196d5b485abSSatish Balay row = idx_i[j]; 197ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 198d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 199d5b485abSSatish Balay ctr[proc]++; 200d5b485abSSatish Balay *ptr[proc] = row; 201d5b485abSSatish Balay ptr[proc]++; 202d6b45a43SBarry Smith } else { /* Update the local table */ 203f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 204d5b485abSSatish Balay } 205d5b485abSSatish Balay } 206d5b485abSSatish Balay /* Update the headers for the current IS */ 207d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 208d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 209d5b485abSSatish Balay outdat_j = outdat[j]; 210d5b485abSSatish Balay k = ++outdat_j[0]; 211d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 212d5b485abSSatish Balay outdat_j[2*k-1] = i; 213d5b485abSSatish Balay } 214d5b485abSSatish Balay } 215d5b485abSSatish Balay isz[i] = isz_i; 216d5b485abSSatish Balay } 217d5b485abSSatish Balay } 218d5b485abSSatish Balay 219d5b485abSSatish Balay /* Now post the sends */ 220b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 221d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 222d5b485abSSatish Balay j = pa[i]; 223b24ad042SBarry Smith ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 224d5b485abSSatish Balay } 225d5b485abSSatish Balay 226d5b485abSSatish Balay /* No longer need the original indices*/ 227d5b485abSSatish Balay for (i=0; i<imax; ++i) { 228d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 229d5b485abSSatish Balay } 23024c049a4SHong Zhang ierr = PetscFree2(idx,n);CHKERRQ(ierr); 231d5b485abSSatish Balay 232d5b485abSSatish Balay for (i=0; i<imax; ++i) { 2336bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 234d5b485abSSatish Balay } 235d5b485abSSatish Balay 236d5b485abSSatish Balay /* Do Local work*/ 237df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 238d5b485abSSatish Balay 239d5b485abSSatish Balay /* Receive messages*/ 240b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr); 2410c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);} 242d5b485abSSatish Balay 243b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 2440c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 245d5b485abSSatish Balay 246d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 24723bdbc58SBarry Smith ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr); 24823bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 249d5b485abSSatish Balay 250b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&xdata);CHKERRQ(ierr); 251b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt),&isz1);CHKERRQ(ierr); 252df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 253c05d87d6SBarry Smith ierr = PetscFree(rbuf[0]);CHKERRQ(ierr); 254606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 255d5b485abSSatish Balay 256d5b485abSSatish Balay /* Send the data back*/ 257d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 258d5b485abSSatish Balay { 259b24ad042SBarry Smith PetscMPIInt *rw1; 260d5b485abSSatish Balay 261b24ad042SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&rw1);CHKERRQ(ierr); 262b24ad042SBarry Smith ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr); 263c7dd2924SSatish Balay 264d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 265d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 266e32f2f54SBarry Smith if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 267d5b485abSSatish Balay rw1[proc] = isz1[i]; 268d5b485abSSatish Balay } 269d5b485abSSatish Balay 270c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 271c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 272d5b485abSSatish Balay 273c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 274c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 27503512dc5SSatish Balay ierr = PetscFree(rw1);CHKERRQ(ierr); 276c7dd2924SSatish Balay } 277c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 278c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 279d5b485abSSatish Balay 280d5b485abSSatish Balay /* Now post the sends */ 281b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 282d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 283d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 284b24ad042SBarry Smith ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 285d5b485abSSatish Balay } 286d5b485abSSatish Balay 287d5b485abSSatish Balay /* receive work done on other processors*/ 288d5b485abSSatish Balay { 289b24ad042SBarry Smith PetscMPIInt idex; 290b24ad042SBarry Smith PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 291f1af5d2fSBarry Smith PetscBT table_i; 292d5b485abSSatish Balay MPI_Status *status2; 293d5b485abSSatish Balay 294169449a0SSatish Balay ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr); 295d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 296999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 297d5b485abSSatish Balay /* Process the message*/ 298999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 299d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 300999d9058SBarry Smith jmax = rbuf2[idex][0]; 301d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 302d5b485abSSatish Balay max = rbuf2_i[2*j]; 303d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 304d5b485abSSatish Balay isz_i = isz[is_no]; 305d5b485abSSatish Balay data_i = data[is_no]; 306d5b485abSSatish Balay table_i = table[is_no]; 307d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 308d5b485abSSatish Balay row = rbuf2_i[ct1]; 309f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 310d5b485abSSatish Balay } 311d5b485abSSatish Balay isz[is_no] = isz_i; 312d5b485abSSatish Balay } 313d5b485abSSatish Balay } 3140c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);} 315606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 316d5b485abSSatish Balay } 317d5b485abSSatish Balay 318d5b485abSSatish Balay for (i=0; i<imax; ++i) { 31970b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 320d5b485abSSatish Balay } 321d5b485abSSatish Balay 322c7dd2924SSatish Balay 323c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 324c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 325c7dd2924SSatish Balay 326606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 327c05d87d6SBarry Smith ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr); 328606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 329606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 330606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 331606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 332606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 3338f9f447aSHong Zhang ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr); 334606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 335606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 336606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 337606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 338606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 3393a40ed3dSBarry Smith PetscFunctionReturn(0); 340d5b485abSSatish Balay } 341d5b485abSSatish Balay 3424a2ae208SSatish Balay #undef __FUNCT__ 3434a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local" 344d5b485abSSatish Balay /* 345df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 346d5b485abSSatish Balay the work on the local processor. 347d5b485abSSatish Balay 348d5b485abSSatish Balay Inputs: 349df36731bSBarry Smith C - MAT_MPIBAIJ; 350d5b485abSSatish Balay imax - total no of index sets processed at a time; 351df36731bSBarry Smith table - an array of char - size = Mbs bits. 352d5b485abSSatish Balay 353d5b485abSSatish Balay Output: 3540e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 355d5b485abSSatish Balay to each index set; 356d5b485abSSatish Balay data - pointer to the solutions 357d5b485abSSatish Balay */ 358b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 359d5b485abSSatish Balay { 360df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 361d5b485abSSatish Balay Mat A = c->A,B = c->B; 362df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 363b24ad042SBarry Smith PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 364b24ad042SBarry Smith PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 365f1af5d2fSBarry Smith PetscBT table_i; 366d5b485abSSatish Balay 3673a40ed3dSBarry Smith PetscFunctionBegin; 368899cda47SBarry Smith rstart = c->rstartbs; 369899cda47SBarry Smith cstart = c->cstartbs; 370d5b485abSSatish Balay ai = a->i; 371df36731bSBarry Smith aj = a->j; 372d5b485abSSatish Balay bi = b->i; 373df36731bSBarry Smith bj = b->j; 374d5b485abSSatish Balay garray = c->garray; 375d5b485abSSatish Balay 376d5b485abSSatish Balay 377d5b485abSSatish Balay for (i=0; i<imax; i++) { 378d5b485abSSatish Balay data_i = data[i]; 379d5b485abSSatish Balay table_i = table[i]; 380d5b485abSSatish Balay isz_i = isz[i]; 381d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 382d5b485abSSatish Balay row = data_i[j] - rstart; 383d5b485abSSatish Balay start = ai[row]; 384d5b485abSSatish Balay end = ai[row+1]; 385d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 386df36731bSBarry Smith val = aj[k] + cstart; 387f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 388d5b485abSSatish Balay } 389d5b485abSSatish Balay start = bi[row]; 390d5b485abSSatish Balay end = bi[row+1]; 391d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 392df36731bSBarry Smith val = garray[bj[k]]; 393f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 394d5b485abSSatish Balay } 395d5b485abSSatish Balay } 396d5b485abSSatish Balay isz[i] = isz_i; 397d5b485abSSatish Balay } 3983a40ed3dSBarry Smith PetscFunctionReturn(0); 399d5b485abSSatish Balay } 4004a2ae208SSatish Balay #undef __FUNCT__ 4014a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive" 402d5b485abSSatish Balay /* 403df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 404d5b485abSSatish Balay and return the output 405d5b485abSSatish Balay 406d5b485abSSatish Balay Input: 407d5b485abSSatish Balay C - the matrix 408d5b485abSSatish Balay nrqr - no of messages being processed. 409d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 410d5b485abSSatish Balay 411d5b485abSSatish Balay Output: 412d5b485abSSatish Balay xdata - array of messages to be sent back 413d5b485abSSatish Balay isz1 - size of each message 414d5b485abSSatish Balay 415a8c7a070SBarry Smith For better efficiency perhaps we should malloc separately each xdata[i], 416d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 4170e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of 418d5b485abSSatish Balay memory is used. 419d5b485abSSatish Balay 420d5b485abSSatish Balay */ 421b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 422d5b485abSSatish Balay { 423df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 424d5b485abSSatish Balay Mat A = c->A,B = c->B; 425df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 4266849ba73SBarry Smith PetscErrorCode ierr; 427b24ad042SBarry Smith PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 428b24ad042SBarry Smith PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 429d2a212eaSBarry Smith PetscInt val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 430b24ad042SBarry Smith PetscInt *rbuf_i,kmax,rbuf_0; 431f1af5d2fSBarry Smith PetscBT xtable; 432d5b485abSSatish Balay 4333a40ed3dSBarry Smith PetscFunctionBegin; 434df36731bSBarry Smith Mbs = c->Mbs; 435899cda47SBarry Smith rstart = c->rstartbs; 436899cda47SBarry Smith cstart = c->cstartbs; 437d5b485abSSatish Balay ai = a->i; 438df36731bSBarry Smith aj = a->j; 439d5b485abSSatish Balay bi = b->i; 440df36731bSBarry Smith bj = b->j; 441d5b485abSSatish Balay garray = c->garray; 442d5b485abSSatish Balay 443d5b485abSSatish Balay 444d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 445d5b485abSSatish Balay rbuf_i = rbuf[i]; 446d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 447d5b485abSSatish Balay ct += rbuf_0; 448d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; } 449d5b485abSSatish Balay } 450d5b485abSSatish Balay 451701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 452701b8814SBarry Smith else max1 = 1; 453d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 454b24ad042SBarry Smith ierr = PetscMalloc(mem_estimate*sizeof(PetscInt),&xdata[0]);CHKERRQ(ierr); 455d5b485abSSatish Balay ++no_malloc; 45653b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr); 457b24ad042SBarry Smith ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr); 458d5b485abSSatish Balay 459d5b485abSSatish Balay ct3 = 0; 460d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 461d5b485abSSatish Balay rbuf_i = rbuf[i]; 462d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 463d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 464d5b485abSSatish Balay ct2 = ct1; 465d5b485abSSatish Balay ct3 += ct1; 466d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 4676831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 468d5b485abSSatish Balay oct2 = ct2; 469d5b485abSSatish Balay kmax = rbuf_i[2*j]; 470d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 471d5b485abSSatish Balay row = rbuf_i[ct1]; 472f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 473d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 474b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 475b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 476b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 477606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 478d5b485abSSatish Balay xdata[0] = tmp; 479d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 480d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 481d5b485abSSatish Balay } 482d5b485abSSatish Balay xdata[i][ct2++] = row; 483d5b485abSSatish Balay ct3++; 484d5b485abSSatish Balay } 485d5b485abSSatish Balay } 486d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 487d5b485abSSatish Balay row = xdata[i][k] - rstart; 488d5b485abSSatish Balay start = ai[row]; 489d5b485abSSatish Balay end = ai[row+1]; 490d5b485abSSatish Balay for (l=start; l<end; l++) { 491df36731bSBarry Smith val = aj[l] + cstart; 492f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 493d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 494b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 495b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 496b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 497606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 498d5b485abSSatish Balay xdata[0] = tmp; 499d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 500d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 501d5b485abSSatish Balay } 502d5b485abSSatish Balay xdata[i][ct2++] = val; 503d5b485abSSatish Balay ct3++; 504d5b485abSSatish Balay } 505d5b485abSSatish Balay } 506d5b485abSSatish Balay start = bi[row]; 507d5b485abSSatish Balay end = bi[row+1]; 508d5b485abSSatish Balay for (l=start; l<end; l++) { 509df36731bSBarry Smith val = garray[bj[l]]; 510f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 511d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 512b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 513b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 514b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 515606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 516d5b485abSSatish Balay xdata[0] = tmp; 517d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 518d5b485abSSatish Balay for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 519d5b485abSSatish Balay } 520d5b485abSSatish Balay xdata[i][ct2++] = val; 521d5b485abSSatish Balay ct3++; 522d5b485abSSatish Balay } 523d5b485abSSatish Balay } 524d5b485abSSatish Balay } 525d5b485abSSatish Balay /* Update the header*/ 526d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 527d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 528d5b485abSSatish Balay } 529d5b485abSSatish Balay xdata[i][0] = rbuf_0; 530d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 531d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 532d5b485abSSatish Balay } 53394bacf5dSBarry Smith ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr); 5341e2582c4SBarry Smith ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr); 5353a40ed3dSBarry Smith PetscFunctionReturn(0); 536d5b485abSSatish Balay } 537d5b485abSSatish Balay 5384a2ae208SSatish Balay #undef __FUNCT__ 5394a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ" 540b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 541d5b485abSSatish Balay { 54236f4e84dSSatish Balay IS *isrow_new,*iscol_new; 543cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 5446849ba73SBarry Smith PetscErrorCode ierr; 5454da72fa9SHong Zhang PetscInt nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs; 5464da72fa9SHong Zhang PetscBool colflag,*allcolumns,*allrows; 547a2feddc5SSatish Balay 5483a40ed3dSBarry Smith PetscFunctionBegin; 54929dcf524SDmitry Karpeev /* Currently, unsorted column indices will result in inverted column indices in the resulting submatrices. */ 55029dcf524SDmitry Karpeev for(i = 0; i < ismax; ++i) { 55129dcf524SDmitry Karpeev PetscBool sorted; 55229dcf524SDmitry Karpeev ierr = ISSorted(iscol[i], &sorted); CHKERRQ(ierr); 553*c7e6e2c7SJed Brown if(!sorted) SETERRQ1(((PetscObject)iscol[i])->comm, PETSC_ERR_SUP, "Column index set %D not sorted", i); 55429dcf524SDmitry Karpeev } 55529dcf524SDmitry Karpeev /* The compression and expansion should be avoided. Doesn't point 55629dcf524SDmitry Karpeev out errors, might change the indices, hence buggey */ 55723bdbc58SBarry Smith ierr = PetscMalloc2(ismax+1,IS,&isrow_new,ismax+1,IS,&iscol_new);CHKERRQ(ierr); 55805d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_new);CHKERRQ(ierr); 55905d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_new);CHKERRQ(ierr); 560cf4f527aSSatish Balay 561307b7a18SHong Zhang /* Check for special case: each processor gets entire matrix columns */ 5624da72fa9SHong Zhang ierr = PetscMalloc2(ismax+1,PetscBool,&allcolumns,ismax+1,PetscBool,&allrows);CHKERRQ(ierr); 563307b7a18SHong Zhang for (i=0; i<ismax; i++) { 564307b7a18SHong Zhang ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 565307b7a18SHong Zhang ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr); 566307b7a18SHong Zhang if (colflag && ncol == C->cmap->N){ 567307b7a18SHong Zhang allcolumns[i] = PETSC_TRUE; 568307b7a18SHong Zhang } else { 569307b7a18SHong Zhang allcolumns[i] = PETSC_FALSE; 570307b7a18SHong Zhang } 5714da72fa9SHong Zhang 5724da72fa9SHong Zhang ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 5734da72fa9SHong Zhang ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr); 5744da72fa9SHong Zhang if (colflag && nrow == C->rmap->N){ 5754da72fa9SHong Zhang allrows[i] = PETSC_TRUE; 5764da72fa9SHong Zhang } else { 5774da72fa9SHong Zhang allrows[i] = PETSC_FALSE; 5784da72fa9SHong Zhang } 579307b7a18SHong Zhang } 580307b7a18SHong Zhang 581cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 582cf4f527aSSatish Balay if (scall != MAT_REUSE_MATRIX) { 58382502324SSatish Balay ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr); 584cf4f527aSSatish Balay } 585cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 586b24ad042SBarry Smith nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 587329f5518SBarry Smith if (!nmax) nmax = 1; 588cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 589cf4f527aSSatish Balay 590653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 5917adad957SLisandro Dalcin ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);CHKERRQ(ierr); 592cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 593cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 594cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 595cf4f527aSSatish Balay else max_no = ismax-pos; 5964da72fa9SHong Zhang ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 597cf4f527aSSatish Balay pos += max_no; 598cf4f527aSSatish Balay } 59936f4e84dSSatish Balay 60036f4e84dSSatish Balay for (i=0; i<ismax; i++) { 6016bf464f9SBarry Smith ierr = ISDestroy(&isrow_new[i]);CHKERRQ(ierr); 6026bf464f9SBarry Smith ierr = ISDestroy(&iscol_new[i]);CHKERRQ(ierr); 60336f4e84dSSatish Balay } 60423bdbc58SBarry Smith ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr); 6054da72fa9SHong Zhang ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr); 6063a40ed3dSBarry Smith PetscFunctionReturn(0); 607a2feddc5SSatish Balay } 608a2feddc5SSatish Balay 609233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 6104a2ae208SSatish Balay #undef __FUNCT__ 6114a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc" 612e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 613233c2ff6SSatish Balay { 614e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 615b9f7ace7SBarry Smith PetscMPIInt fproc = PetscMPIIntCast( (PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5))); 616233c2ff6SSatish Balay 617233c2ff6SSatish Balay PetscFunctionBegin; 61823ce1328SBarry Smith if (fproc > size) fproc = size; 619e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 620e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 621233c2ff6SSatish Balay else fproc++; 622233c2ff6SSatish Balay } 623e005ede5SBarry Smith *rank = fproc; 624233c2ff6SSatish Balay PetscFunctionReturn(0); 625233c2ff6SSatish Balay } 626233c2ff6SSatish Balay #endif 627233c2ff6SSatish Balay 628a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 629b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 6304a2ae208SSatish Balay #undef __FUNCT__ 6314a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local" 6324da72fa9SHong Zhang PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 633a2feddc5SSatish Balay { 634df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 635cf4f527aSSatish Balay Mat A = c->A; 636df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 6375d0c19d7SBarry Smith const PetscInt **irow,**icol,*irow_i; 6385d0c19d7SBarry Smith PetscInt *nrow,*ncol,*w3,*w4,start; 6396849ba73SBarry Smith PetscErrorCode ierr; 6409405f653SBarry Smith PetscMPIInt size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc; 6419405f653SBarry Smith PetscInt **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row; 642b24ad042SBarry Smith PetscInt nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol; 643b24ad042SBarry Smith PetscInt **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 6445d0c19d7SBarry Smith PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 645052f0c41SBarry Smith PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 646d0f46423SBarry Smith PetscInt bs=C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 647899cda47SBarry Smith PetscInt cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 648899cda47SBarry Smith PetscInt *bmap = c->garray,ctmp,rstart=c->rstartbs; 6497a868f3eSHong Zhang MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3; 6507a868f3eSHong Zhang MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3; 651d5b485abSSatish Balay MPI_Comm comm; 652ace3abfcSBarry Smith PetscBool flag; 653b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1; 6547a868f3eSHong Zhang PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 6557a868f3eSHong Zhang /* variables below are used for the matrix numerical values - case of !ijonly */ 6567a868f3eSHong Zhang MPI_Request *r_waits4,*s_waits4; 6577a868f3eSHong Zhang MPI_Status *r_status4,*s_status4; 658f3e30f2cSMatthew G Knepley MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a = PETSC_NULL,*sbuf_aa_i,*vworkA = PETSC_NULL,*vworkB = PETSC_NULL; 6597a868f3eSHong Zhang MatScalar *a_a=a->a,*b_a=b->a; 660c7dd2924SSatish Balay 66180d608c0SSatish Balay #if defined (PETSC_USE_CTABLE) 662b24ad042SBarry Smith PetscInt tt; 6639066e3d5SHong Zhang PetscTable *rmap,*cmap,rmap_i,cmap_i=PETSC_NULL; 664233c2ff6SSatish Balay #else 6659066e3d5SHong Zhang PetscInt **cmap,*cmap_i=PETSC_NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 666233c2ff6SSatish Balay #endif 667d5b485abSSatish Balay 6683a40ed3dSBarry Smith PetscFunctionBegin; 6697adad957SLisandro Dalcin comm = ((PetscObject)C)->comm; 6707adad957SLisandro Dalcin tag0 = ((PetscObject)C)->tag; 671d5b485abSSatish Balay size = c->size; 672d5b485abSSatish Balay rank = c->rank; 673d5b485abSSatish Balay 67434e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 67534e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 67634e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 67734e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 67834e6ae18SSatish Balay 679052f0c41SBarry Smith #if defined(PETSC_USE_CTABLE) 68023bdbc58SBarry Smith ierr = PetscMalloc4(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol);CHKERRQ(ierr); 681052f0c41SBarry Smith #else 68223bdbc58SBarry Smith ierr = PetscMalloc5(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol,Mbs+1,PetscInt,&rtable);CHKERRQ(ierr); 683d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 684d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 68582a7e548SBarry Smith jmax = C->rmap->range[i+1]/bs; 686d5b485abSSatish Balay for (; j<jmax; j++) { 687d5b485abSSatish Balay rtable[j] = i; 688d5b485abSSatish Balay } 689d5b485abSSatish Balay } 690233c2ff6SSatish Balay #endif 691233c2ff6SSatish Balay 692233c2ff6SSatish Balay for (i=0; i<ismax; i++) { 6934da72fa9SHong Zhang if (allrows[i]){ 6944da72fa9SHong Zhang irow[i] = PETSC_NULL; 6954da72fa9SHong Zhang nrow[i] = C->rmap->N/bs; 6964da72fa9SHong Zhang } else { 697233c2ff6SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 698233c2ff6SSatish Balay ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 6994da72fa9SHong Zhang } 7004da72fa9SHong Zhang 701307b7a18SHong Zhang if (allcolumns[i]){ 702307b7a18SHong Zhang icol[i] = PETSC_NULL; 703307b7a18SHong Zhang ncol[i] = C->cmap->N/bs; 704307b7a18SHong Zhang } else { 705307b7a18SHong Zhang ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 706233c2ff6SSatish Balay ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 707233c2ff6SSatish Balay } 708307b7a18SHong Zhang } 709d5b485abSSatish Balay 710d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg,and buffer space 711d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 71223bdbc58SBarry Smith ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr); 71323bdbc58SBarry Smith ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 71423bdbc58SBarry Smith ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 71523bdbc58SBarry Smith ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr); 716d5b485abSSatish Balay for (i=0; i<ismax; i++) { 717b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 718d5b485abSSatish Balay jmax = nrow[i]; 719d5b485abSSatish Balay irow_i = irow[i]; 720d5b485abSSatish Balay for (j=0; j<jmax; j++) { 7214da72fa9SHong Zhang if (allrows[i]) { 7224da72fa9SHong Zhang row = j; 7234da72fa9SHong Zhang } else { 724d5b485abSSatish Balay row = irow_i[j]; 7254da72fa9SHong Zhang } 726233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 727899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 728233c2ff6SSatish Balay #else 729d5b485abSSatish Balay proc = rtable[row]; 730233c2ff6SSatish Balay #endif 731d5b485abSSatish Balay w4[proc]++; 732d5b485abSSatish Balay } 733d5b485abSSatish Balay for (j=0; j<size; j++) { 734d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 735d5b485abSSatish Balay } 736d5b485abSSatish Balay } 737d5b485abSSatish Balay 738d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 739e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 740d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 741d5b485abSSatish Balay w3[rank] = 0; 742d5b485abSSatish Balay for (i=0; i<size; i++) { 743d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 744d5b485abSSatish Balay } 745b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/ 746d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 747d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 748d5b485abSSatish Balay } 749d5b485abSSatish Balay 750d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 751d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 752d5b485abSSatish Balay j = pa[i]; 753d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 754d5b485abSSatish Balay msz += w1[j]; 755d5b485abSSatish Balay } 756d5b485abSSatish Balay 757c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 758a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 759c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 760d5b485abSSatish Balay 761c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 762c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 763c7dd2924SSatish Balay 764c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 765c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 766d5b485abSSatish Balay 767d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 76823bdbc58SBarry Smith ierr = PetscMalloc4(size,PetscInt*,&sbuf1,size,PetscInt*,&ptr,2*msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr); 76923bdbc58SBarry Smith ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 77023bdbc58SBarry Smith ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 771d5b485abSSatish Balay { 772b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 773d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 774d5b485abSSatish Balay j = pa[i]; 775d5b485abSSatish Balay iptr += ict; 776d5b485abSSatish Balay sbuf1[j] = iptr; 777d5b485abSSatish Balay ict = w1[j]; 778d5b485abSSatish Balay } 779d5b485abSSatish Balay } 780d5b485abSSatish Balay 781d5b485abSSatish Balay /* Form the outgoing messages */ 782d5b485abSSatish Balay /* Initialise the header space */ 783d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 784d5b485abSSatish Balay j = pa[i]; 785d5b485abSSatish Balay sbuf1[j][0] = 0; 786b24ad042SBarry Smith ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 787d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 788d5b485abSSatish Balay } 789d5b485abSSatish Balay 790d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 791d5b485abSSatish Balay for (i=0; i<ismax; i++) { 792b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 793d5b485abSSatish Balay irow_i = irow[i]; 794d5b485abSSatish Balay jmax = nrow[i]; 795d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 7964da72fa9SHong Zhang if (allrows[i]){ 7974da72fa9SHong Zhang row = j; 7984da72fa9SHong Zhang } else { 799d5b485abSSatish Balay row = irow_i[j]; 8004da72fa9SHong Zhang } 801233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 802899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 803233c2ff6SSatish Balay #else 804d5b485abSSatish Balay proc = rtable[row]; 805233c2ff6SSatish Balay #endif 806d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 807d5b485abSSatish Balay ctr[proc]++; 808d5b485abSSatish Balay *ptr[proc] = row; 809d5b485abSSatish Balay ptr[proc]++; 810d5b485abSSatish Balay } 811d5b485abSSatish Balay } 812d5b485abSSatish Balay /* Update the headers for the current IS */ 813d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 81406ef35abSSatish Balay if ((ctr_j = ctr[j])) { 815d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 816d5b485abSSatish Balay k = ++sbuf1_j[0]; 817d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 818d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 819d5b485abSSatish Balay } 820d5b485abSSatish Balay } 821d5b485abSSatish Balay } 822d5b485abSSatish Balay 823d5b485abSSatish Balay /* Now post the sends */ 824b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 825d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 826d5b485abSSatish Balay j = pa[i]; 827b24ad042SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 828d5b485abSSatish Balay } 829d5b485abSSatish Balay 830d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 831b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr); 832b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr); 833d5b485abSSatish Balay rbuf2[0] = tmp + msz; 834d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 835d5b485abSSatish Balay j = pa[i]; 836d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 837d5b485abSSatish Balay } 838d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 839d5b485abSSatish Balay j = pa[i]; 840b24ad042SBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 841d5b485abSSatish Balay } 842d5b485abSSatish Balay 843d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 844d5b485abSSatish Balay 845d5b485abSSatish Balay /* Receive messages*/ 846b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 847b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr); 84823bdbc58SBarry Smith ierr = PetscMalloc3(nrqr+1,PetscInt*,&sbuf2,nrqr,PetscInt,&req_size,nrqr,PetscInt,&req_source);CHKERRQ(ierr); 849d5b485abSSatish Balay { 850df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 851b24ad042SBarry Smith PetscInt *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 852d5b485abSSatish Balay 853d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 854999d9058SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 855999d9058SBarry Smith req_size[idex] = 0; 856999d9058SBarry Smith rbuf1_i = rbuf1[idex]; 857d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 858b24ad042SBarry Smith ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 859b24ad042SBarry Smith ierr = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr); 860999d9058SBarry Smith sbuf2_i = sbuf2[idex]; 861d5b485abSSatish Balay for (j=start; j<end; j++) { 862d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 863d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 864d5b485abSSatish Balay sbuf2_i[j] = ncols; 865999d9058SBarry Smith req_size[idex] += ncols; 866d5b485abSSatish Balay } 867999d9058SBarry Smith req_source[idex] = r_status1[i].MPI_SOURCE; 868d5b485abSSatish Balay /* form the header */ 869999d9058SBarry Smith sbuf2_i[0] = req_size[idex]; 870d5b485abSSatish Balay for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; } 871b24ad042SBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 872d5b485abSSatish Balay } 873d5b485abSSatish Balay } 874606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 875606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 876d5b485abSSatish Balay 877d5b485abSSatish Balay /* recv buffer sizes */ 878d5b485abSSatish Balay /* Receive messages*/ 879b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr); 880b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr); 881b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr); 8827a868f3eSHong Zhang if (!ijonly){ 8837a868f3eSHong Zhang ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr); 8847a868f3eSHong Zhang ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr); 8857a868f3eSHong Zhang } 886d5b485abSSatish Balay 887d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 888999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 889b24ad042SBarry Smith ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr); 890b24ad042SBarry Smith ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 8917a868f3eSHong Zhang if (!ijonly){ 8927a868f3eSHong Zhang ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr); 893b24ad042SBarry Smith ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 894d5b485abSSatish Balay } 8957a868f3eSHong Zhang } 896606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 897606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 898d5b485abSSatish Balay 899d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 900b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr); 901b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr); 902d5b485abSSatish Balay 9030c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 9040c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 905606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 906606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 907606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 908606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 909d5b485abSSatish Balay 910d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 911b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr); 912d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 913b24ad042SBarry Smith ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr); 914d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 915d5b485abSSatish Balay 916b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr); 917d5b485abSSatish Balay { 918d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 919d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 920d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 921d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 922d5b485abSSatish Balay ct2 = 0; 923d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 924d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 925d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 926e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 927d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 928d5b485abSSatish Balay ncols = nzA + nzB; 929d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 930d5b485abSSatish Balay 931d5b485abSSatish Balay /* load the column indices for this row into cols*/ 932d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 933d5b485abSSatish Balay for (l=0; l<nzB; l++) { 934d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 935d5b485abSSatish Balay else break; 936d5b485abSSatish Balay } 937d5b485abSSatish Balay imark = l; 938d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 939d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 940d5b485abSSatish Balay ct2 += ncols; 941d5b485abSSatish Balay } 942d5b485abSSatish Balay } 943b24ad042SBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 944d5b485abSSatish Balay } 945d5b485abSSatish Balay } 946b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr); 947b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr); 948d5b485abSSatish Balay 949d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 9507a868f3eSHong Zhang if (!ijonly){ 95182502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr); 952d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 95382502324SSatish Balay ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr); 954a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 955d5b485abSSatish Balay 956b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr); 957d5b485abSSatish Balay { 958d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 959d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 960d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 961d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 962d5b485abSSatish Balay ct2 = 0; 963d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 964d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 965d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 966e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 967d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 968d5b485abSSatish Balay ncols = nzA + nzB; 969d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 970a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 971d5b485abSSatish Balay 972d5b485abSSatish Balay /* load the column values for this row into vals*/ 9735b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 974d5b485abSSatish Balay for (l=0; l<nzB; l++) { 9753a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 9763eda8832SBarry Smith ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9773a40ed3dSBarry Smith } 978d5b485abSSatish Balay else break; 979d5b485abSSatish Balay } 980d5b485abSSatish Balay imark = l; 9813a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 9823eda8832SBarry Smith ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9833a40ed3dSBarry Smith } 9843a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 9853eda8832SBarry Smith ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9863a40ed3dSBarry Smith } 987d5b485abSSatish Balay ct2 += ncols; 988d5b485abSSatish Balay } 989d5b485abSSatish Balay } 9903eda8832SBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 991d5b485abSSatish Balay } 992d5b485abSSatish Balay } 993b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr); 994b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr); 9957a868f3eSHong Zhang } 996533163c2SBarry Smith ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 997606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 998d5b485abSSatish Balay 999d5b485abSSatish Balay /* Form the matrix */ 1000307b7a18SHong Zhang /* create col map: global col of C -> local col of submatrices */ 1001d5b485abSSatish Balay { 10025d0c19d7SBarry Smith const PetscInt *icol_i; 1003233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 10048fa52d88SHong Zhang ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&cmap);CHKERRQ(ierr); 1005ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 1006307b7a18SHong Zhang if (!allcolumns[i]){ 10078fa52d88SHong Zhang ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 1008307b7a18SHong Zhang jmax = ncol[i]; 1009307b7a18SHong Zhang icol_i = icol[i]; 10108fa52d88SHong Zhang cmap_i = cmap[i]; 1011307b7a18SHong Zhang for (j=0; j<jmax; j++) { 10123861aac3SJed Brown ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1013307b7a18SHong Zhang } 1014307b7a18SHong Zhang } else { 10158fa52d88SHong Zhang cmap[i] = PETSC_NULL; 1016307b7a18SHong Zhang } 1017233c2ff6SSatish Balay } 1018233c2ff6SSatish Balay #else 101923bdbc58SBarry Smith ierr = PetscMalloc(ismax*sizeof(PetscInt*),&cmap);CHKERRQ(ierr); 1020d5b485abSSatish Balay for (i=0; i<ismax; i++) { 10218fa52d88SHong Zhang if (!allcolumns[i]){ 10228fa52d88SHong Zhang ierr = PetscMalloc(c->Nbs*sizeof(PetscInt),&cmap[i]);CHKERRQ(ierr); 10238fa52d88SHong Zhang ierr = PetscMemzero(cmap[i],c->Nbs*sizeof(PetscInt));CHKERRQ(ierr); 1024d5b485abSSatish Balay jmax = ncol[i]; 1025d5b485abSSatish Balay icol_i = icol[i]; 1026d5b485abSSatish Balay cmap_i = cmap[i]; 1027d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1028d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1029d5b485abSSatish Balay } 10308fa52d88SHong Zhang } else { /* allcolumns[i] */ 10318fa52d88SHong Zhang cmap[i] = PETSC_NULL; 10328fa52d88SHong Zhang } 1033d5b485abSSatish Balay } 1034307b7a18SHong Zhang #endif 1035d5b485abSSatish Balay } 1036d5b485abSSatish Balay 1037d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 1038d5b485abSSatish Balay for (i=0,j=0; i<ismax; i++) { j += nrow[i]; } 1039052f0c41SBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr); 1040b24ad042SBarry Smith lens[0] = (PetscInt*)(lens + ismax); 1041b24ad042SBarry Smith ierr = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr); 1042d5b485abSSatish Balay for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; } 1043d5b485abSSatish Balay 1044d5b485abSSatish Balay /* Update lens from local data */ 1045d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1046d5b485abSSatish Balay jmax = nrow[i]; 10478fa52d88SHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 1048d5b485abSSatish Balay irow_i = irow[i]; 1049d5b485abSSatish Balay lens_i = lens[i]; 1050d5b485abSSatish Balay for (j=0; j<jmax; j++) { 10514da72fa9SHong Zhang if (allrows[i]){ 10524da72fa9SHong Zhang row = j; 10534da72fa9SHong Zhang } else { 1054d5b485abSSatish Balay row = irow_i[j]; 10554da72fa9SHong Zhang } 1056233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1057899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1058233c2ff6SSatish Balay #else 1059d5b485abSSatish Balay proc = rtable[row]; 1060233c2ff6SSatish Balay #endif 1061d5b485abSSatish Balay if (proc == rank) { 106236f4e84dSSatish Balay /* Get indices from matA and then from matB */ 106336f4e84dSSatish Balay row = row - rstart; 106436f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 106536f4e84dSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1066307b7a18SHong Zhang if (!allcolumns[i]){ 1067233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1068233c2ff6SSatish Balay for (k=0; k<nzA; k++) { 10698fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1070233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1071233c2ff6SSatish Balay } 1072233c2ff6SSatish Balay for (k=0; k<nzB; k++) { 10738fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1074233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1075233c2ff6SSatish Balay } 1076307b7a18SHong Zhang 1077233c2ff6SSatish Balay #else 1078ca161407SBarry Smith for (k=0; k<nzA; k++) { 107936f4e84dSSatish Balay if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; } 1080ca161407SBarry Smith } 1081ca161407SBarry Smith for (k=0; k<nzB; k++) { 108236f4e84dSSatish Balay if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; } 1083d5b485abSSatish Balay } 1084233c2ff6SSatish Balay #endif 1085307b7a18SHong Zhang } else {/* allcolumns */ 1086307b7a18SHong Zhang lens_i[j] = nzA + nzB; 1087307b7a18SHong Zhang } 1088d5b485abSSatish Balay } 1089a2feddc5SSatish Balay } 1090ca161407SBarry Smith } 1091233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1092d5b485abSSatish Balay /* Create row map*/ 10938fa52d88SHong Zhang ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rmap);CHKERRQ(ierr); 1094ff0794f7SSatish Balay for (i=0; i<ismax; i++){ 10958fa52d88SHong Zhang ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1096233c2ff6SSatish Balay } 1097233c2ff6SSatish Balay #else 1098233c2ff6SSatish Balay /* Create row map*/ 1099052f0c41SBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr); 1100b24ad042SBarry Smith rmap[0] = (PetscInt*)(rmap + ismax); 1101b24ad042SBarry Smith ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr); 1102233c2ff6SSatish Balay for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;} 1103233c2ff6SSatish Balay #endif 1104d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1105d5b485abSSatish Balay irow_i = irow[i]; 1106d5b485abSSatish Balay jmax = nrow[i]; 1107233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 11088fa52d88SHong Zhang rmap_i = rmap[i]; 1109233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 11104da72fa9SHong Zhang if (allrows[i]){ 11113861aac3SJed Brown ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 11124da72fa9SHong Zhang } else { 11133861aac3SJed Brown ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1114233c2ff6SSatish Balay } 11154da72fa9SHong Zhang } 1116233c2ff6SSatish Balay #else 1117233c2ff6SSatish Balay rmap_i = rmap[i]; 1118d5b485abSSatish Balay for (j=0; j<jmax; j++) { 11194da72fa9SHong Zhang if (allrows[i]){ 11204da72fa9SHong Zhang rmap_i[j] = j; 11214da72fa9SHong Zhang } else { 1122d5b485abSSatish Balay rmap_i[irow_i[j]] = j; 1123d5b485abSSatish Balay } 11244da72fa9SHong Zhang } 1125233c2ff6SSatish Balay #endif 1126d5b485abSSatish Balay } 1127d5b485abSSatish Balay 1128d5b485abSSatish Balay /* Update lens from offproc data */ 1129d5b485abSSatish Balay { 1130b24ad042SBarry Smith PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1131b24ad042SBarry Smith PetscMPIInt ii; 1132d5b485abSSatish Balay 1133d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1134b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr); 1135b24ad042SBarry Smith idex = pa[ii]; 1136999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1137d5b485abSSatish Balay jmax = sbuf1_i[0]; 1138d5b485abSSatish Balay ct1 = 2*jmax+1; 1139d5b485abSSatish Balay ct2 = 0; 1140b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1141b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 1142d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1143d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1144d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1145d5b485abSSatish Balay lens_i = lens[is_no]; 11468fa52d88SHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1147d5b485abSSatish Balay rmap_i = rmap[is_no]; 1148d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1149233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 11508fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1151233c2ff6SSatish Balay row--; 1152e32f2f54SBarry Smith if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); } 1153233c2ff6SSatish Balay #else 1154d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1155233c2ff6SSatish Balay #endif 1156d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1157d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1158307b7a18SHong Zhang if (!allcolumns[is_no]){ 1159233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 11608fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 1161233c2ff6SSatish Balay if (tt) { 1162233c2ff6SSatish Balay lens_i[row]++; 1163233c2ff6SSatish Balay } 1164233c2ff6SSatish Balay #else 1165d5b485abSSatish Balay if (cmap_i[rbuf3_i[ct2]]) { 1166d5b485abSSatish Balay lens_i[row]++; 1167d5b485abSSatish Balay } 1168233c2ff6SSatish Balay #endif 1169307b7a18SHong Zhang } else { /* allcolumns */ 1170307b7a18SHong Zhang lens_i[row]++; 1171307b7a18SHong Zhang } 1172d5b485abSSatish Balay } 1173d5b485abSSatish Balay } 1174d5b485abSSatish Balay } 1175d5b485abSSatish Balay } 1176d5b485abSSatish Balay } 1177606d414cSSatish Balay ierr = PetscFree(r_status3);CHKERRQ(ierr); 1178606d414cSSatish Balay ierr = PetscFree(r_waits3);CHKERRQ(ierr); 11790c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 1180606d414cSSatish Balay ierr = PetscFree(s_status3);CHKERRQ(ierr); 1181606d414cSSatish Balay ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1182d5b485abSSatish Balay 1183d5b485abSSatish Balay /* Create the submatrices */ 1184d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 11857a868f3eSHong Zhang if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet"); 1186d5b485abSSatish Balay /* 1187d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1188d5b485abSSatish Balay */ 1189d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1190df36731bSBarry Smith mat = (Mat_SeqBAIJ *)(submats[i]->data); 1191e7e72b3dSBarry Smith if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 1192b24ad042SBarry Smith ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 1193e7e72b3dSBarry Smith if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 1194d5b485abSSatish Balay /* Initial matrix as if empty */ 1195b24ad042SBarry Smith ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr); 1196d5f3da31SBarry Smith submats[i]->factortype = C->factortype; 1197d5b485abSSatish Balay } 1198ca161407SBarry Smith } else { 11997a868f3eSHong Zhang PetscInt bs_tmp; 12007a868f3eSHong Zhang if (ijonly){ 12017a868f3eSHong Zhang bs_tmp = 1; 12027a868f3eSHong Zhang } else { 12037a868f3eSHong Zhang bs_tmp = bs; 12047a868f3eSHong Zhang } 1205d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1206f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 12077a868f3eSHong Zhang ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr); 12087adad957SLisandro Dalcin ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 12097a868f3eSHong Zhang ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 12107a868f3eSHong Zhang ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 1211d5b485abSSatish Balay } 1212d5b485abSSatish Balay } 1213d5b485abSSatish Balay 1214d5b485abSSatish Balay /* Assemble the matrices */ 1215d5b485abSSatish Balay /* First assemble the local rows */ 1216d5b485abSSatish Balay { 1217b24ad042SBarry Smith PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 1218f3e30f2cSMatthew G Knepley MatScalar *imat_a = PETSC_NULL; 1219d5b485abSSatish Balay 1220d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1221df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[i]->data; 1222d5b485abSSatish Balay imat_ilen = mat->ilen; 1223d5b485abSSatish Balay imat_j = mat->j; 1224d5b485abSSatish Balay imat_i = mat->i; 12257a868f3eSHong Zhang if (!ijonly) imat_a = mat->a; 12268fa52d88SHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 1227d5b485abSSatish Balay rmap_i = rmap[i]; 1228d5b485abSSatish Balay irow_i = irow[i]; 1229d5b485abSSatish Balay jmax = nrow[i]; 1230d5b485abSSatish Balay for (j=0; j<jmax; j++) { 12314da72fa9SHong Zhang if (allrows[i]){ 12324da72fa9SHong Zhang row = j; 12334da72fa9SHong Zhang } else { 1234d5b485abSSatish Balay row = irow_i[j]; 12354da72fa9SHong Zhang } 1236233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1237899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1238233c2ff6SSatish Balay #else 1239d5b485abSSatish Balay proc = rtable[row]; 1240233c2ff6SSatish Balay #endif 1241d5b485abSSatish Balay if (proc == rank) { 124236f4e84dSSatish Balay row = row - rstart; 124336f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 124436f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 124536f4e84dSSatish Balay cworkA = a_j + a_i[row]; 124636f4e84dSSatish Balay cworkB = b_j + b_i[row]; 12477a868f3eSHong Zhang if (!ijonly){ 124836f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 124936f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 12507a868f3eSHong Zhang } 1251233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 12528fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1253233c2ff6SSatish Balay row--; 1254e32f2f54SBarry Smith if (row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); } 1255233c2ff6SSatish Balay #else 125636f4e84dSSatish Balay row = rmap_i[row + rstart]; 1257233c2ff6SSatish Balay #endif 1258df36731bSBarry Smith mat_i = imat_i[row]; 12597a868f3eSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1260d5b485abSSatish Balay mat_j = imat_j + mat_i; 126136f4e84dSSatish Balay ilen_row = imat_ilen[row]; 126236f4e84dSSatish Balay 126336f4e84dSSatish Balay /* load the column indices for this row into cols*/ 1264307b7a18SHong Zhang if (!allcolumns[i]){ 126536f4e84dSSatish Balay for (l=0; l<nzB; l++) { 126636f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 1267233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 12688fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1269233c2ff6SSatish Balay if (tcol) { 1270233c2ff6SSatish Balay #else 127136f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1272233c2ff6SSatish Balay #endif 1273df36731bSBarry Smith *mat_j++ = tcol - 1; 12743eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1275549d3d68SSatish Balay mat_a += bs2; 1276d5b485abSSatish Balay ilen_row++; 1277d5b485abSSatish Balay } 1278ca161407SBarry Smith } else break; 127936f4e84dSSatish Balay } 128036f4e84dSSatish Balay imark = l; 128136f4e84dSSatish Balay for (l=0; l<nzA; l++) { 1282233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 12838fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1284233c2ff6SSatish Balay if (tcol) { 1285233c2ff6SSatish Balay #else 128636f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 1287233c2ff6SSatish Balay #endif 128836f4e84dSSatish Balay *mat_j++ = tcol - 1; 12897a868f3eSHong Zhang if (!ijonly){ 12903eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1291549d3d68SSatish Balay mat_a += bs2; 12927a868f3eSHong Zhang } 129336f4e84dSSatish Balay ilen_row++; 129436f4e84dSSatish Balay } 129536f4e84dSSatish Balay } 129636f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 1297233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 12988fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1299233c2ff6SSatish Balay if (tcol) { 1300233c2ff6SSatish Balay #else 130136f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1302233c2ff6SSatish Balay #endif 130336f4e84dSSatish Balay *mat_j++ = tcol - 1; 13047a868f3eSHong Zhang if (!ijonly){ 13053eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1306549d3d68SSatish Balay mat_a += bs2; 13077a868f3eSHong Zhang } 130836f4e84dSSatish Balay ilen_row++; 130936f4e84dSSatish Balay } 131036f4e84dSSatish Balay } 1311307b7a18SHong Zhang } else { /* allcolumns */ 1312307b7a18SHong Zhang for (l=0; l<nzB; l++) { 1313307b7a18SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1314307b7a18SHong Zhang *mat_j++ = ctmp; 1315307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1316307b7a18SHong Zhang mat_a += bs2; 1317307b7a18SHong Zhang ilen_row++; 1318307b7a18SHong Zhang } else break; 1319307b7a18SHong Zhang } 1320307b7a18SHong Zhang imark = l; 1321307b7a18SHong Zhang for (l=0; l<nzA; l++) { 1322307b7a18SHong Zhang *mat_j++ = cstart+cworkA[l]; 1323307b7a18SHong Zhang if (!ijonly){ 1324307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1325307b7a18SHong Zhang mat_a += bs2; 1326307b7a18SHong Zhang } 1327307b7a18SHong Zhang ilen_row++; 1328307b7a18SHong Zhang } 1329307b7a18SHong Zhang for (l=imark; l<nzB; l++) { 1330307b7a18SHong Zhang *mat_j++ = bmap[cworkB[l]]; 1331307b7a18SHong Zhang if (!ijonly){ 1332307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1333307b7a18SHong Zhang mat_a += bs2; 1334307b7a18SHong Zhang } 1335307b7a18SHong Zhang ilen_row++; 1336307b7a18SHong Zhang } 1337307b7a18SHong Zhang } 1338d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1339d5b485abSSatish Balay } 1340d5b485abSSatish Balay } 1341d5b485abSSatish Balay } 1342d5b485abSSatish Balay } 1343d5b485abSSatish Balay 1344d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1345d5b485abSSatish Balay { 1346b24ad042SBarry Smith PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1347b24ad042SBarry Smith PetscInt *imat_j,*imat_i; 1348f3e30f2cSMatthew G Knepley MatScalar *imat_a = PETSC_NULL,*rbuf4_i = PETSC_NULL; 1349b24ad042SBarry Smith PetscMPIInt ii; 1350d5b485abSSatish Balay 1351d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 13527a868f3eSHong Zhang if (ijonly){ 13537a868f3eSHong Zhang ii = tmp2; 13547a868f3eSHong Zhang } else { 1355b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr); 13567a868f3eSHong Zhang } 1357b24ad042SBarry Smith idex = pa[ii]; 1358999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1359d5b485abSSatish Balay jmax = sbuf1_i[0]; 1360d5b485abSSatish Balay ct1 = 2*jmax + 1; 1361d5b485abSSatish Balay ct2 = 0; 1362b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1363b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 13647a868f3eSHong Zhang if (!ijonly) rbuf4_i = rbuf4[ii]; 1365d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1366d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 13678fa52d88SHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1368d5b485abSSatish Balay rmap_i = rmap[is_no]; 1369df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[is_no]->data; 1370d5b485abSSatish Balay imat_ilen = mat->ilen; 1371d5b485abSSatish Balay imat_j = mat->j; 1372d5b485abSSatish Balay imat_i = mat->i; 13737a868f3eSHong Zhang if (!ijonly) imat_a = mat->a; 1374d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1375d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1376d5b485abSSatish Balay row = sbuf1_i[ct1]; 1377233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 13788fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 1379233c2ff6SSatish Balay row--; 1380e32f2f54SBarry Smith if(row < 0) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); } 1381233c2ff6SSatish Balay #else 1382d5b485abSSatish Balay row = rmap_i[row]; 1383233c2ff6SSatish Balay #endif 1384d5b485abSSatish Balay ilen = imat_ilen[row]; 1385df36731bSBarry Smith mat_i = imat_i[row]; 13867a868f3eSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1387d5b485abSSatish Balay mat_j = imat_j + mat_i; 1388d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1389307b7a18SHong Zhang 1390307b7a18SHong Zhang if (!allcolumns[is_no]){ 1391d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1392233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 13938fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1394233c2ff6SSatish Balay if (tcol) { 1395233c2ff6SSatish Balay #else 1396d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1397233c2ff6SSatish Balay #endif 1398df36731bSBarry Smith *mat_j++ = tcol - 1; 13997a868f3eSHong Zhang if (!ijonly){ 14003eda8832SBarry Smith ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1401549d3d68SSatish Balay mat_a += bs2; 14027a868f3eSHong Zhang } 1403d5b485abSSatish Balay ilen++; 1404d5b485abSSatish Balay } 1405d5b485abSSatish Balay } 1406307b7a18SHong Zhang } else { /* allcolumns */ 1407307b7a18SHong Zhang for (l=0; l<max2; l++,ct2++) { 1408307b7a18SHong Zhang *mat_j++ = rbuf3_i[ct2]; 1409307b7a18SHong Zhang if (!ijonly){ 1410307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1411307b7a18SHong Zhang mat_a += bs2; 1412307b7a18SHong Zhang } 1413307b7a18SHong Zhang ilen++; 1414307b7a18SHong Zhang } 1415307b7a18SHong Zhang } 1416d5b485abSSatish Balay imat_ilen[row] = ilen; 1417d5b485abSSatish Balay } 1418d5b485abSSatish Balay } 1419d5b485abSSatish Balay } 1420d5b485abSSatish Balay } 14217a868f3eSHong Zhang if (!ijonly){ 1422606d414cSSatish Balay ierr = PetscFree(r_status4);CHKERRQ(ierr); 1423606d414cSSatish Balay ierr = PetscFree(r_waits4);CHKERRQ(ierr); 14240c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 1425606d414cSSatish Balay ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1426606d414cSSatish Balay ierr = PetscFree(s_status4);CHKERRQ(ierr); 14277a868f3eSHong Zhang } 1428d5b485abSSatish Balay 1429d5b485abSSatish Balay /* Restore the indices */ 1430d5b485abSSatish Balay for (i=0; i<ismax; i++) { 14314da72fa9SHong Zhang if (!allrows[i]){ 1432d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 14334da72fa9SHong Zhang } 1434307b7a18SHong Zhang if (!allcolumns[i]){ 1435d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1436d5b485abSSatish Balay } 1437307b7a18SHong Zhang } 1438d5b485abSSatish Balay 1439d5b485abSSatish Balay /* Destroy allocated memory */ 144023bdbc58SBarry Smith #if defined(PETSC_USE_CTABLE) 144123bdbc58SBarry Smith ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr); 144223bdbc58SBarry Smith #else 144323bdbc58SBarry Smith ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr); 144423bdbc58SBarry Smith #endif 144523bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 1446606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 1447d5b485abSSatish Balay 144823bdbc58SBarry Smith ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 1449606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 1450606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 1451d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1452606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1453d5b485abSSatish Balay } 1454d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1455606d414cSSatish Balay ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1456d5b485abSSatish Balay } 145723bdbc58SBarry Smith ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr); 1458606d414cSSatish Balay ierr = PetscFree(rbuf3);CHKERRQ(ierr); 1459606d414cSSatish Balay ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1460606d414cSSatish Balay ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 14617a868f3eSHong Zhang if (!ijonly) { 14627a868f3eSHong Zhang for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);} 14637a868f3eSHong Zhang ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1464606d414cSSatish Balay ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1465606d414cSSatish Balay ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 14667a868f3eSHong Zhang } 1467d5b485abSSatish Balay 1468233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1469ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 14708fa52d88SHong Zhang ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr); 1471233c2ff6SSatish Balay } 1472233c2ff6SSatish Balay #endif 14738fa52d88SHong Zhang ierr = PetscFree(rmap);CHKERRQ(ierr); 14748fa52d88SHong Zhang 14758fa52d88SHong Zhang for (i=0; i<ismax; i++){ 14768fa52d88SHong Zhang if (!allcolumns[i]){ 14778fa52d88SHong Zhang #if defined (PETSC_USE_CTABLE) 14788fa52d88SHong Zhang ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr); 14798fa52d88SHong Zhang #else 14808fa52d88SHong Zhang ierr = PetscFree(cmap[i]);CHKERRQ(ierr); 14818fa52d88SHong Zhang #endif 14828fa52d88SHong Zhang } 14838fa52d88SHong Zhang } 14848fa52d88SHong Zhang ierr = PetscFree(cmap);CHKERRQ(ierr); 1485606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1486d5b485abSSatish Balay 1487d5b485abSSatish Balay for (i=0; i<ismax; i++) { 148836f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 148936f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1490d5b485abSSatish Balay } 14917a868f3eSHong Zhang 14927a868f3eSHong Zhang c->ijonly = PETSC_FALSE; /* set back to the default */ 14933a40ed3dSBarry Smith PetscFunctionReturn(0); 1494d5b485abSSatish Balay } 1495ca161407SBarry Smith 1496