1d5b485abSSatish Balay /* 2d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 3d5b485abSSatish Balay and to find submatrices that were shared across processors. 4d5b485abSSatish Balay */ 570f55243SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h" 60a835dfdSSatish Balay #include "petscbt.h" 7d5b485abSSatish Balay 8*b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat,PetscInt,IS *); 9*b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char **,PetscInt*,PetscInt**); 10*b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt **,PetscInt**,PetscInt*); 11*b24ad042SBarry Smith EXTERN PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 12*b24ad042SBarry Smith EXTERN PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 13d5b485abSSatish Balay 144a2ae208SSatish Balay #undef __FUNCT__ 154a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ" 16*b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov) 17d5b485abSSatish Balay { 18d9489beaSHong Zhang Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 196849ba73SBarry Smith PetscErrorCode ierr; 20*b24ad042SBarry Smith PetscInt i,N=C->N, bs=c->bs; 2136f4e84dSSatish Balay IS *is_new; 2236f4e84dSSatish Balay 233a40ed3dSBarry Smith PetscFunctionBegin; 2482502324SSatish Balay ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr); 25df36731bSBarry Smith /* Convert the indices into block format */ 2627f478b1SHong Zhang ierr = ISCompressIndicesGeneral(N,bs,imax,is,is_new);CHKERRQ(ierr); 2729bbc08cSBarry Smith if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");} 28d5b485abSSatish Balay for (i=0; i<ov; ++i) { 2936f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 30d5b485abSSatish Balay } 31ca161407SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);} 3227f478b1SHong Zhang ierr = ISExpandIndicesGeneral(N,bs,imax,is_new,is);CHKERRQ(ierr); 33ca161407SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);} 34606d414cSSatish Balay ierr = PetscFree(is_new);CHKERRQ(ierr); 353a40ed3dSBarry Smith PetscFunctionReturn(0); 36d5b485abSSatish Balay } 37d5b485abSSatish Balay 38d5b485abSSatish Balay /* 39d5b485abSSatish Balay Sample message format: 40d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 410e9b0e7eSHong Zhang to index sets is[1], is[5] 42d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 43d5b485abSSatish Balay ----------- 44d5b485abSSatish Balay mesg [1] = 1 => is[1] 45d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 46d5b485abSSatish Balay ----------- 47d5b485abSSatish Balay mesg [5] = 5 => is[5] 48d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 49d5b485abSSatish Balay ----------- 50d5b485abSSatish Balay mesg [7] 510e9b0e7eSHong Zhang mesg [n] data(is[1]) 52d5b485abSSatish Balay ----------- 53d5b485abSSatish Balay mesg[n+1] 54d5b485abSSatish Balay mesg[m] data(is[5]) 55d5b485abSSatish Balay ----------- 56d5b485abSSatish Balay 57d5b485abSSatish Balay Notes: 58d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 59d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 60d5b485abSSatish Balay */ 614a2ae208SSatish Balay #undef __FUNCT__ 624a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once" 63*b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[]) 64d5b485abSSatish Balay { 65df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 66*b24ad042SBarry Smith PetscInt **idx,*n,*w3,*w4,*rtable,**data,len,*idx_i; 676849ba73SBarry Smith PetscErrorCode ierr; 68*b24ad042SBarry Smith PetscMPIInt size,rank,tag1,tag2,*w2,*w1,nrqr; 69*b24ad042SBarry Smith PetscInt Mbs,i,j,k,**rbuf,row,proc,nrqs,msz,**outdat,**ptr; 70*b24ad042SBarry Smith PetscInt *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2; 71*b24ad042SBarry 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; 76d5b485abSSatish Balay 773a40ed3dSBarry Smith PetscFunctionBegin; 78d5b485abSSatish Balay comm = C->comm; 79d5b485abSSatish Balay size = c->size; 80d5b485abSSatish Balay rank = c->rank; 81df36731bSBarry Smith Mbs = c->Mbs; 82d5b485abSSatish Balay 83c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 84c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 85c7dd2924SSatish Balay 86*b24ad042SBarry Smith len = (imax+1)*sizeof(PetscInt*)+ (imax + Mbs)*sizeof(PetscInt); 8782502324SSatish Balay ierr = PetscMalloc(len,&idx);CHKERRQ(ierr); 88*b24ad042SBarry Smith n = (PetscInt*)(idx + imax); 89d5b485abSSatish Balay rtable = n + imax; 90d5b485abSSatish Balay 91d5b485abSSatish Balay for (i=0; i<imax; i++) { 92d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 93b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 94d5b485abSSatish Balay } 95d5b485abSSatish Balay 96d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 97d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 98d5b485abSSatish Balay len = c->rowners[i+1]; 99d5b485abSSatish Balay for (; j<len; j++) { 100d5b485abSSatish Balay rtable[j] = i; 101d5b485abSSatish Balay } 102d5b485abSSatish Balay } 103d5b485abSSatish Balay 104d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 105d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 106*b24ad042SBarry Smith ierr = PetscMalloc(size*2*sizeof(PetscInt)+size*2*sizeof(PetscMPIInt),&w1);CHKERRQ(ierr);/* mesg size */ 107d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 108*b24ad042SBarry Smith w3 = (PetscInt*) (w2 + size); /* no of IS that needs to be sent to proc i */ 109d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 110*b24ad042SBarry Smith ierr = PetscMemzero(w1,size*sizeof(PetscInt)+2*size*sizeof(PetscMPIInt));CHKERRQ(ierr); /* initialise work vector*/ 111d5b485abSSatish Balay for (i=0; i<imax; i++) { 112*b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 113d5b485abSSatish Balay idx_i = idx[i]; 114d5b485abSSatish Balay len = n[i]; 115d5b485abSSatish Balay for (j=0; j<len; j++) { 116d5b485abSSatish Balay row = idx_i[j]; 1176b41c64dSBarry Smith if (row < 0) { 118e005ede5SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 1196b41c64dSBarry Smith } 120d5b485abSSatish Balay proc = rtable[row]; 121d5b485abSSatish Balay w4[proc]++; 122d5b485abSSatish Balay } 123d5b485abSSatish Balay for (j=0; j<size; j++){ 124d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 125d5b485abSSatish Balay } 126d5b485abSSatish Balay } 127d5b485abSSatish Balay 128d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 129d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1300e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 131d5b485abSSatish Balay w3[rank] = 0; 132d5b485abSSatish Balay for (i=0; i<size; i++) { 133d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 134d5b485abSSatish Balay } 135d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 136*b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); 137d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 138d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 139d5b485abSSatish Balay } 140d5b485abSSatish Balay 141d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 142d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 143d5b485abSSatish Balay j = pa[i]; 144d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 145d5b485abSSatish Balay msz += w1[j]; 146d5b485abSSatish Balay } 147d5b485abSSatish Balay 148c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 149a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 150c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 151d5b485abSSatish Balay 152c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 153c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 154d5b485abSSatish Balay 155d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 156*b24ad042SBarry Smith len = 2*size*sizeof(PetscInt*) + (size+msz)*sizeof(PetscInt); 15782502324SSatish Balay ierr = PetscMalloc(len,&outdat);CHKERRQ(ierr); 158d5b485abSSatish Balay ptr = outdat + size; /* Pointers to the data in outgoing buffers */ 159*b24ad042SBarry Smith ierr = PetscMemzero(outdat,2*size*sizeof(PetscInt*));CHKERRQ(ierr); 160*b24ad042SBarry Smith tmp = (PetscInt*)(outdat + 2*size); 161d5b485abSSatish Balay ctr = tmp + msz; 162d5b485abSSatish Balay 163d5b485abSSatish Balay { 164*b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 165d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 166d5b485abSSatish Balay j = pa[i]; 167d5b485abSSatish Balay iptr += ict; 168d5b485abSSatish Balay outdat[j] = iptr; 169d5b485abSSatish Balay ict = w1[j]; 170d5b485abSSatish Balay } 171d5b485abSSatish Balay } 172d5b485abSSatish Balay 173d5b485abSSatish Balay /* Form the outgoing messages */ 174d5b485abSSatish Balay /*plug in the headers*/ 175d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 176d5b485abSSatish Balay j = pa[i]; 177d5b485abSSatish Balay outdat[j][0] = 0; 178*b24ad042SBarry Smith ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 179d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 180d5b485abSSatish Balay } 181d5b485abSSatish Balay 182d5b485abSSatish Balay /* Memory for doing local proc's work*/ 183d5b485abSSatish Balay { 184*b24ad042SBarry Smith PetscInt *d_p; 185d5b485abSSatish Balay char *t_p; 186d5b485abSSatish Balay 187*b24ad042SBarry Smith len = (imax)*(sizeof(PetscBT) + sizeof(PetscInt*)+ sizeof(PetscInt)) + 188*b24ad042SBarry Smith (Mbs)*imax*sizeof(PetscInt) + (Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char) + 1; 18982502324SSatish Balay ierr = PetscMalloc(len,&table);CHKERRQ(ierr); 190549d3d68SSatish Balay ierr = PetscMemzero(table,len);CHKERRQ(ierr); 191*b24ad042SBarry Smith data = (PetscInt **)(table + imax); 192*b24ad042SBarry Smith isz = (PetscInt *)(data + imax); 193*b24ad042SBarry Smith d_p = (PetscInt *)(isz + imax); 194bbd702dbSSatish Balay t_p = (char *)(d_p + Mbs*imax); 195bbd702dbSSatish Balay for (i=0; i<imax; i++) { 196b6410449SSatish Balay table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 197bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 198d5b485abSSatish Balay } 199d5b485abSSatish Balay } 200d5b485abSSatish Balay 201d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 202d5b485abSSatish Balay { 203*b24ad042SBarry Smith PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 204f1af5d2fSBarry Smith PetscBT table_i; 205d5b485abSSatish Balay 206d5b485abSSatish Balay for (i=0; i<imax; i++) { 207*b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 208d5b485abSSatish Balay n_i = n[i]; 209d5b485abSSatish Balay table_i = table[i]; 210d5b485abSSatish Balay idx_i = idx[i]; 211d5b485abSSatish Balay data_i = data[i]; 212d5b485abSSatish Balay isz_i = isz[i]; 213d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 214d5b485abSSatish Balay row = idx_i[j]; 215d5b485abSSatish Balay proc = rtable[row]; 216d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 217d5b485abSSatish Balay ctr[proc]++; 218d5b485abSSatish Balay *ptr[proc] = row; 219d5b485abSSatish Balay ptr[proc]++; 220d5b485abSSatish Balay } 221d5b485abSSatish Balay else { /* Update the local table */ 222f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 223d5b485abSSatish Balay } 224d5b485abSSatish Balay } 225d5b485abSSatish Balay /* Update the headers for the current IS */ 226d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 227d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 228d5b485abSSatish Balay outdat_j = outdat[j]; 229d5b485abSSatish Balay k = ++outdat_j[0]; 230d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 231d5b485abSSatish Balay outdat_j[2*k-1] = i; 232d5b485abSSatish Balay } 233d5b485abSSatish Balay } 234d5b485abSSatish Balay isz[i] = isz_i; 235d5b485abSSatish Balay } 236d5b485abSSatish Balay } 237d5b485abSSatish Balay 238d5b485abSSatish Balay /* Now post the sends */ 239b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 240d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 241d5b485abSSatish Balay j = pa[i]; 242*b24ad042SBarry Smith ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 243d5b485abSSatish Balay } 244d5b485abSSatish Balay 245d5b485abSSatish Balay /* No longer need the original indices*/ 246d5b485abSSatish Balay for (i=0; i<imax; ++i) { 247d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 248d5b485abSSatish Balay } 249606d414cSSatish Balay ierr = PetscFree(idx);CHKERRQ(ierr); 250d5b485abSSatish Balay 251d5b485abSSatish Balay for (i=0; i<imax; ++i) { 252d5b485abSSatish Balay ierr = ISDestroy(is[i]);CHKERRQ(ierr); 253d5b485abSSatish Balay } 254d5b485abSSatish Balay 255d5b485abSSatish Balay /* Do Local work*/ 256df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 257d5b485abSSatish Balay 258d5b485abSSatish Balay /* Receive messages*/ 259b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr); 260c7dd2924SSatish Balay ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr); 261d5b485abSSatish Balay 262b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 263ca161407SBarry Smith ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr); 264d5b485abSSatish Balay 265d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 266606d414cSSatish Balay ierr = PetscFree(outdat);CHKERRQ(ierr); 267606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 268d5b485abSSatish Balay 269*b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&xdata);CHKERRQ(ierr); 270*b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt),&isz1);CHKERRQ(ierr); 271df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 272606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 273d5b485abSSatish Balay 274d5b485abSSatish Balay /* Send the data back*/ 275d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 276d5b485abSSatish Balay { 277*b24ad042SBarry Smith PetscMPIInt *rw1; 278d5b485abSSatish Balay 279*b24ad042SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&rw1);CHKERRQ(ierr); 280*b24ad042SBarry Smith ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr); 281c7dd2924SSatish Balay 282d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 283d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 284e005ede5SBarry Smith if (proc != onodes1[i]) SETERRQ(PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 285d5b485abSSatish Balay rw1[proc] = isz1[i]; 286d5b485abSSatish Balay } 287d5b485abSSatish Balay 288c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 289c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 290d5b485abSSatish Balay 291c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 292c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 293c7dd2924SSatish Balay PetscFree(rw1); 294c7dd2924SSatish Balay } 295c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 296c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 297d5b485abSSatish Balay 298d5b485abSSatish Balay /* Now post the sends */ 299b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 300d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 301d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 302*b24ad042SBarry Smith ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 303d5b485abSSatish Balay } 304d5b485abSSatish Balay 305d5b485abSSatish Balay /* receive work done on other processors*/ 306d5b485abSSatish Balay { 307*b24ad042SBarry Smith PetscMPIInt idex; 308*b24ad042SBarry Smith PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 309f1af5d2fSBarry Smith PetscBT table_i; 310d5b485abSSatish Balay MPI_Status *status2; 311d5b485abSSatish Balay 312169449a0SSatish Balay ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr); 313d5b485abSSatish Balay 314d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 315999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 316d5b485abSSatish Balay /* Process the message*/ 317999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 318d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 319999d9058SBarry Smith jmax = rbuf2[idex][0]; 320d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 321d5b485abSSatish Balay max = rbuf2_i[2*j]; 322d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 323d5b485abSSatish Balay isz_i = isz[is_no]; 324d5b485abSSatish Balay data_i = data[is_no]; 325d5b485abSSatish Balay table_i = table[is_no]; 326d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 327d5b485abSSatish Balay row = rbuf2_i[ct1]; 328f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 329d5b485abSSatish Balay } 330d5b485abSSatish Balay isz[is_no] = isz_i; 331d5b485abSSatish Balay } 332d5b485abSSatish Balay } 333ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr); 334606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 335d5b485abSSatish Balay } 336d5b485abSSatish Balay 337d5b485abSSatish Balay for (i=0; i<imax; ++i) { 338029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr); 339d5b485abSSatish Balay } 340d5b485abSSatish Balay 341c7dd2924SSatish Balay 342c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 343c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 344c7dd2924SSatish Balay 345606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 346606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 347606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 348606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 349606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 350606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 351606d414cSSatish Balay ierr = PetscFree(table);CHKERRQ(ierr); 352606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 353606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 354606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 355606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 356606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 3573a40ed3dSBarry Smith PetscFunctionReturn(0); 358d5b485abSSatish Balay } 359d5b485abSSatish Balay 3604a2ae208SSatish Balay #undef __FUNCT__ 3614a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local" 362d5b485abSSatish Balay /* 363df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 364d5b485abSSatish Balay the work on the local processor. 365d5b485abSSatish Balay 366d5b485abSSatish Balay Inputs: 367df36731bSBarry Smith C - MAT_MPIBAIJ; 368d5b485abSSatish Balay imax - total no of index sets processed at a time; 369df36731bSBarry Smith table - an array of char - size = Mbs bits. 370d5b485abSSatish Balay 371d5b485abSSatish Balay Output: 3720e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 373d5b485abSSatish Balay to each index set; 374d5b485abSSatish Balay data - pointer to the solutions 375d5b485abSSatish Balay */ 376*b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 377d5b485abSSatish Balay { 378df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 379d5b485abSSatish Balay Mat A = c->A,B = c->B; 380df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 381*b24ad042SBarry Smith PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 382*b24ad042SBarry Smith PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 383f1af5d2fSBarry Smith PetscBT table_i; 384d5b485abSSatish Balay 3853a40ed3dSBarry Smith PetscFunctionBegin; 386d5b485abSSatish Balay rstart = c->rstart; 387d5b485abSSatish Balay cstart = c->cstart; 388d5b485abSSatish Balay ai = a->i; 389df36731bSBarry Smith aj = a->j; 390d5b485abSSatish Balay bi = b->i; 391df36731bSBarry Smith bj = b->j; 392d5b485abSSatish Balay garray = c->garray; 393d5b485abSSatish Balay 394d5b485abSSatish Balay 395d5b485abSSatish Balay for (i=0; i<imax; i++) { 396d5b485abSSatish Balay data_i = data[i]; 397d5b485abSSatish Balay table_i = table[i]; 398d5b485abSSatish Balay isz_i = isz[i]; 399d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 400d5b485abSSatish Balay row = data_i[j] - rstart; 401d5b485abSSatish Balay start = ai[row]; 402d5b485abSSatish Balay end = ai[row+1]; 403d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 404df36731bSBarry Smith val = aj[k] + cstart; 405f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 406d5b485abSSatish Balay } 407d5b485abSSatish Balay start = bi[row]; 408d5b485abSSatish Balay end = bi[row+1]; 409d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 410df36731bSBarry Smith val = garray[bj[k]]; 411f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 412d5b485abSSatish Balay } 413d5b485abSSatish Balay } 414d5b485abSSatish Balay isz[i] = isz_i; 415d5b485abSSatish Balay } 4163a40ed3dSBarry Smith PetscFunctionReturn(0); 417d5b485abSSatish Balay } 4184a2ae208SSatish Balay #undef __FUNCT__ 4194a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive" 420d5b485abSSatish Balay /* 421df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 422d5b485abSSatish Balay and return the output 423d5b485abSSatish Balay 424d5b485abSSatish Balay Input: 425d5b485abSSatish Balay C - the matrix 426d5b485abSSatish Balay nrqr - no of messages being processed. 427d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 428d5b485abSSatish Balay 429d5b485abSSatish Balay Output: 430d5b485abSSatish Balay xdata - array of messages to be sent back 431d5b485abSSatish Balay isz1 - size of each message 432d5b485abSSatish Balay 433d5b485abSSatish Balay For better efficiency perhaps we should malloc seperately each xdata[i], 434d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 4350e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of 436d5b485abSSatish Balay memory is used. 437d5b485abSSatish Balay 438d5b485abSSatish Balay */ 439*b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 440d5b485abSSatish Balay { 441df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 442d5b485abSSatish Balay Mat A = c->A,B = c->B; 443df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 4446849ba73SBarry Smith PetscErrorCode ierr; 445*b24ad042SBarry Smith PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 446*b24ad042SBarry Smith PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 447*b24ad042SBarry Smith PetscInt val,max1,max2,rank,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 448*b24ad042SBarry Smith PetscInt *rbuf_i,kmax,rbuf_0; 449f1af5d2fSBarry Smith PetscBT xtable; 450d5b485abSSatish Balay 4513a40ed3dSBarry Smith PetscFunctionBegin; 452d5b485abSSatish Balay rank = c->rank; 453df36731bSBarry Smith Mbs = c->Mbs; 454d5b485abSSatish Balay rstart = c->rstart; 455d5b485abSSatish Balay cstart = c->cstart; 456d5b485abSSatish Balay ai = a->i; 457df36731bSBarry Smith aj = a->j; 458d5b485abSSatish Balay bi = b->i; 459df36731bSBarry Smith bj = b->j; 460d5b485abSSatish Balay garray = c->garray; 461d5b485abSSatish Balay 462d5b485abSSatish Balay 463d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 464d5b485abSSatish Balay rbuf_i = rbuf[i]; 465d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 466d5b485abSSatish Balay ct += rbuf_0; 467d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; } 468d5b485abSSatish Balay } 469d5b485abSSatish Balay 470701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 471701b8814SBarry Smith else max1 = 1; 472d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 473*b24ad042SBarry Smith ierr = PetscMalloc(mem_estimate*sizeof(PetscInt),&xdata[0]);CHKERRQ(ierr); 474d5b485abSSatish Balay ++no_malloc; 4756831982aSBarry Smith ierr = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr); 476*b24ad042SBarry Smith ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr); 477d5b485abSSatish Balay 478d5b485abSSatish Balay ct3 = 0; 479d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 480d5b485abSSatish Balay rbuf_i = rbuf[i]; 481d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 482d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 483d5b485abSSatish Balay ct2 = ct1; 484d5b485abSSatish Balay ct3 += ct1; 485d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 4866831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 487d5b485abSSatish Balay oct2 = ct2; 488d5b485abSSatish Balay kmax = rbuf_i[2*j]; 489d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 490d5b485abSSatish Balay row = rbuf_i[ct1]; 491f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 492d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 493*b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 494*b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 495*b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 496606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 497d5b485abSSatish Balay xdata[0] = tmp; 498d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 499d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 500d5b485abSSatish Balay } 501d5b485abSSatish Balay xdata[i][ct2++] = row; 502d5b485abSSatish Balay ct3++; 503d5b485abSSatish Balay } 504d5b485abSSatish Balay } 505d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 506d5b485abSSatish Balay row = xdata[i][k] - rstart; 507d5b485abSSatish Balay start = ai[row]; 508d5b485abSSatish Balay end = ai[row+1]; 509d5b485abSSatish Balay for (l=start; l<end; l++) { 510df36731bSBarry Smith val = aj[l] + cstart; 511f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 512d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 513*b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 514*b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 515*b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 516606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 517d5b485abSSatish Balay xdata[0] = tmp; 518d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 519d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 520d5b485abSSatish Balay } 521d5b485abSSatish Balay xdata[i][ct2++] = val; 522d5b485abSSatish Balay ct3++; 523d5b485abSSatish Balay } 524d5b485abSSatish Balay } 525d5b485abSSatish Balay start = bi[row]; 526d5b485abSSatish Balay end = bi[row+1]; 527d5b485abSSatish Balay for (l=start; l<end; l++) { 528df36731bSBarry Smith val = garray[bj[l]]; 529f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 530d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 531*b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 532*b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 533*b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 534606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 535d5b485abSSatish Balay xdata[0] = tmp; 536d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 537d5b485abSSatish Balay for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 538d5b485abSSatish Balay } 539d5b485abSSatish Balay xdata[i][ct2++] = val; 540d5b485abSSatish Balay ct3++; 541d5b485abSSatish Balay } 542d5b485abSSatish Balay } 543d5b485abSSatish Balay } 544d5b485abSSatish Balay /* Update the header*/ 545d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 546d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 547d5b485abSSatish Balay } 548d5b485abSSatish Balay xdata[i][0] = rbuf_0; 549d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 550d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 551d5b485abSSatish Balay } 5526831982aSBarry Smith ierr = PetscBTDestroy(xtable);CHKERRQ(ierr); 553b0a32e0cSBarry Smith PetscLogInfo(0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %d bytes, required %d, no of mallocs = %d\n",rank,mem_estimate,ct3,no_malloc); 5543a40ed3dSBarry Smith PetscFunctionReturn(0); 555d5b485abSSatish Balay } 556d5b485abSSatish Balay 557*b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *); 558a2feddc5SSatish Balay 5594a2ae208SSatish Balay #undef __FUNCT__ 5604a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ" 561*b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 562d5b485abSSatish Balay { 56336f4e84dSSatish Balay IS *isrow_new,*iscol_new; 564cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 5656849ba73SBarry Smith PetscErrorCode ierr; 566*b24ad042SBarry Smith PetscInt nmax,nstages_local,nstages,i,pos,max_no,N=C->N,bs=c->bs; 567a2feddc5SSatish Balay 5683a40ed3dSBarry Smith PetscFunctionBegin; 5693a40ed3dSBarry Smith /* The compression and expansion should be avoided. Does'nt point 570a2feddc5SSatish Balay out errors might change the indices hence buggey */ 571a2feddc5SSatish Balay 572b0a32e0cSBarry Smith ierr = PetscMalloc(2*(ismax+1)*sizeof(IS),&isrow_new);CHKERRQ(ierr); 57336f4e84dSSatish Balay iscol_new = isrow_new + ismax; 57427f478b1SHong Zhang ierr = ISCompressIndicesSorted(N,bs,ismax,isrow,isrow_new);CHKERRQ(ierr); 57527f478b1SHong Zhang ierr = ISCompressIndicesSorted(N,bs,ismax,iscol,iscol_new);CHKERRQ(ierr); 576cf4f527aSSatish Balay 577cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 578cf4f527aSSatish Balay if (scall != MAT_REUSE_MATRIX) { 57982502324SSatish Balay ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr); 580cf4f527aSSatish Balay } 581cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 582*b24ad042SBarry Smith nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 583329f5518SBarry Smith if (!nmax) nmax = 1; 584cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 585cf4f527aSSatish Balay 586653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 587*b24ad042SBarry Smith ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,C->comm);CHKERRQ(ierr); 588cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 589cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 590cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 591cf4f527aSSatish Balay else max_no = ismax-pos; 592cf4f527aSSatish Balay ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr); 593cf4f527aSSatish Balay pos += max_no; 594cf4f527aSSatish Balay } 59536f4e84dSSatish Balay 59636f4e84dSSatish Balay for (i=0; i<ismax; i++) { 597ca161407SBarry Smith ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr); 598ca161407SBarry Smith ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr); 59936f4e84dSSatish Balay } 600606d414cSSatish Balay ierr = PetscFree(isrow_new);CHKERRQ(ierr); 6013a40ed3dSBarry Smith PetscFunctionReturn(0); 602a2feddc5SSatish Balay } 603a2feddc5SSatish Balay 604233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 6054a2ae208SSatish Balay #undef __FUNCT__ 6064a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc" 607e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 608233c2ff6SSatish Balay { 609e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 610*b24ad042SBarry Smith PetscMPIInt fproc = (PetscMPIInt) ((float)row * (float)size / (float)nGlobalNd + 0.5); 611233c2ff6SSatish Balay 612233c2ff6SSatish Balay PetscFunctionBegin; 61323ce1328SBarry Smith if (fproc > size) fproc = size; 614e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 615e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 616233c2ff6SSatish Balay else fproc++; 617233c2ff6SSatish Balay } 618e005ede5SBarry Smith *rank = fproc; 619233c2ff6SSatish Balay PetscFunctionReturn(0); 620233c2ff6SSatish Balay } 621233c2ff6SSatish Balay #endif 622233c2ff6SSatish Balay 623a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 624*b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 6254a2ae208SSatish Balay #undef __FUNCT__ 6264a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local" 627*b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats) 628a2feddc5SSatish Balay { 629df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 630cf4f527aSSatish Balay Mat A = c->A; 631df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 632*b24ad042SBarry Smith PetscInt **irow,**icol,*nrow,*ncol,*w3,*w4,start; 6336849ba73SBarry Smith PetscErrorCode ierr; 634*b24ad042SBarry Smith PetscMPIInt size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end; 635*b24ad042SBarry Smith PetscInt **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row,proc; 636*b24ad042SBarry Smith PetscInt nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol; 637*b24ad042SBarry Smith PetscInt **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 638*b24ad042SBarry Smith PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i; 639*b24ad042SBarry Smith PetscInt len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 640*b24ad042SBarry Smith PetscInt bs=c->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 641*b24ad042SBarry Smith PetscInt cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 642*b24ad042SBarry Smith PetscInt *bmap = c->garray,ctmp,rstart=c->rstart; 643d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 644d5b485abSSatish Balay MPI_Request *r_waits4,*s_waits3,*s_waits4; 645d5b485abSSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 646d5b485abSSatish Balay MPI_Status *r_status3,*r_status4,*s_status4; 647d5b485abSSatish Balay MPI_Comm comm; 6483eda8832SBarry Smith MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB; 6493eda8832SBarry Smith MatScalar *a_a=a->a,*b_a=b->a; 6500f5bd95cSBarry Smith PetscTruth flag; 651*b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1; 652c7dd2924SSatish Balay 65380d608c0SSatish Balay #if defined (PETSC_USE_CTABLE) 654*b24ad042SBarry Smith PetscInt tt; 655233c2ff6SSatish Balay PetscTable *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1; 656233c2ff6SSatish Balay #else 657*b24ad042SBarry Smith PetscInt **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 658233c2ff6SSatish Balay #endif 659d5b485abSSatish Balay 6603a40ed3dSBarry Smith PetscFunctionBegin; 661d5b485abSSatish Balay comm = C->comm; 66234e6ae18SSatish Balay tag0 = C->tag; 663d5b485abSSatish Balay size = c->size; 664d5b485abSSatish Balay rank = c->rank; 665d5b485abSSatish Balay 66634e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 66734e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 66834e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 66934e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 67034e6ae18SSatish Balay 671d5b485abSSatish Balay /* Check if the col indices are sorted */ 672d5b485abSSatish Balay for (i=0; i<ismax; i++) { 673888f2ed8SSatish Balay ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr); 67429bbc08cSBarry Smith if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted"); 675d5b485abSSatish Balay } 676d5b485abSSatish Balay 677*b24ad042SBarry Smith len = (2*ismax+1)*(sizeof(PetscInt*)+ sizeof(PetscInt)); 678233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE) 679*b24ad042SBarry Smith len += (Mbs+1)*sizeof(PetscInt); 680233c2ff6SSatish Balay #endif 68182502324SSatish Balay ierr = PetscMalloc(len,&irow);CHKERRQ(ierr); 682d5b485abSSatish Balay icol = irow + ismax; 683*b24ad042SBarry Smith nrow = (PetscInt*)(icol + ismax); 684d5b485abSSatish Balay ncol = nrow + ismax; 685233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE) 686d5b485abSSatish Balay rtable = ncol + ismax; 687d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 688d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 689d5b485abSSatish Balay jmax = c->rowners[i+1]; 690d5b485abSSatish Balay for (; j<jmax; j++) { 691d5b485abSSatish Balay rtable[j] = i; 692d5b485abSSatish Balay } 693d5b485abSSatish Balay } 694233c2ff6SSatish Balay #endif 695233c2ff6SSatish Balay 696233c2ff6SSatish Balay for (i=0; i<ismax; i++) { 697233c2ff6SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 698233c2ff6SSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 699233c2ff6SSatish Balay ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 700233c2ff6SSatish Balay ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 701233c2ff6SSatish Balay } 702d5b485abSSatish Balay 703d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg,and buffer space 704d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 705*b24ad042SBarry Smith ierr = PetscMalloc(size*2*sizeof(PetscInt) + size*2*sizeof(PetscMPIInt),&w1);CHKERRQ(ierr); /* mesg size */ 706d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 707*b24ad042SBarry Smith w3 = (PetscInt*)(w2 + size); /* no of IS that needs to be sent to proc i */ 708d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 709*b24ad042SBarry Smith ierr = PetscMemzero(w1,size*sizeof(PetscInt)+2*size*sizeof(PetscMPIInt));CHKERRQ(ierr); /* initialise work vector*/ 710d5b485abSSatish Balay for (i=0; i<ismax; i++) { 711*b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 712d5b485abSSatish Balay jmax = nrow[i]; 713d5b485abSSatish Balay irow_i = irow[i]; 714d5b485abSSatish Balay for (j=0; j<jmax; j++) { 715d5b485abSSatish Balay row = irow_i[j]; 716233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 717233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 718233c2ff6SSatish Balay #else 719d5b485abSSatish Balay proc = rtable[row]; 720233c2ff6SSatish Balay #endif 721d5b485abSSatish Balay w4[proc]++; 722d5b485abSSatish Balay } 723d5b485abSSatish Balay for (j=0; j<size; j++) { 724d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 725d5b485abSSatish Balay } 726d5b485abSSatish Balay } 727d5b485abSSatish Balay 728d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 729e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 730d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 731d5b485abSSatish Balay w3[rank] = 0; 732d5b485abSSatish Balay for (i=0; i<size; i++) { 733d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 734d5b485abSSatish Balay } 735*b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/ 736d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 737d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 738d5b485abSSatish Balay } 739d5b485abSSatish Balay 740d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 741d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 742d5b485abSSatish Balay j = pa[i]; 743d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 744d5b485abSSatish Balay msz += w1[j]; 745d5b485abSSatish Balay } 746d5b485abSSatish Balay 747c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 748a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 749c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 750d5b485abSSatish Balay 751c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 752c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 753c7dd2924SSatish Balay 754c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 755c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 756d5b485abSSatish Balay 757d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 758*b24ad042SBarry Smith len = 2*size*sizeof(PetscInt*) + 2*msz*sizeof(PetscInt) + size*sizeof(PetscInt); 75982502324SSatish Balay ierr = PetscMalloc(len,&sbuf1);CHKERRQ(ierr); 760d5b485abSSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 761*b24ad042SBarry Smith ierr = PetscMemzero(sbuf1,2*size*sizeof(PetscInt*));CHKERRQ(ierr); 762d5b485abSSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 763*b24ad042SBarry Smith tmp = (PetscInt*)(ptr + size); 764d5b485abSSatish Balay ctr = tmp + 2*msz; 765d5b485abSSatish Balay 766d5b485abSSatish Balay { 767*b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 768d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 769d5b485abSSatish Balay j = pa[i]; 770d5b485abSSatish Balay iptr += ict; 771d5b485abSSatish Balay sbuf1[j] = iptr; 772d5b485abSSatish Balay ict = w1[j]; 773d5b485abSSatish Balay } 774d5b485abSSatish Balay } 775d5b485abSSatish Balay 776d5b485abSSatish Balay /* Form the outgoing messages */ 777d5b485abSSatish Balay /* Initialise the header space */ 778d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 779d5b485abSSatish Balay j = pa[i]; 780d5b485abSSatish Balay sbuf1[j][0] = 0; 781*b24ad042SBarry Smith ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 782d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 783d5b485abSSatish Balay } 784d5b485abSSatish Balay 785d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 786d5b485abSSatish Balay for (i=0; i<ismax; i++) { 787*b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 788d5b485abSSatish Balay irow_i = irow[i]; 789d5b485abSSatish Balay jmax = nrow[i]; 790d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 791d5b485abSSatish Balay row = irow_i[j]; 792233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 793233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 794233c2ff6SSatish Balay #else 795d5b485abSSatish Balay proc = rtable[row]; 796233c2ff6SSatish Balay #endif 797d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 798d5b485abSSatish Balay ctr[proc]++; 799d5b485abSSatish Balay *ptr[proc] = row; 800d5b485abSSatish Balay ptr[proc]++; 801d5b485abSSatish Balay } 802d5b485abSSatish Balay } 803d5b485abSSatish Balay /* Update the headers for the current IS */ 804d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 80506ef35abSSatish Balay if ((ctr_j = ctr[j])) { 806d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 807d5b485abSSatish Balay k = ++sbuf1_j[0]; 808d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 809d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 810d5b485abSSatish Balay } 811d5b485abSSatish Balay } 812d5b485abSSatish Balay } 813d5b485abSSatish Balay 814d5b485abSSatish Balay /* Now post the sends */ 815b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 816d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 817d5b485abSSatish Balay j = pa[i]; 818*b24ad042SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 819d5b485abSSatish Balay } 820d5b485abSSatish Balay 821d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 822b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr); 823*b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr); 824d5b485abSSatish Balay rbuf2[0] = tmp + msz; 825d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 826d5b485abSSatish Balay j = pa[i]; 827d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 828d5b485abSSatish Balay } 829d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 830d5b485abSSatish Balay j = pa[i]; 831*b24ad042SBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 832d5b485abSSatish Balay } 833d5b485abSSatish Balay 834d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 835d5b485abSSatish Balay 836d5b485abSSatish Balay /* Receive messages*/ 837b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 838b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr); 839*b24ad042SBarry Smith len = 2*nrqr*sizeof(PetscInt) + (nrqr+1)*sizeof(PetscInt*); 84082502324SSatish Balay ierr = PetscMalloc(len,&sbuf2);CHKERRQ(ierr); 841*b24ad042SBarry Smith req_size = (PetscInt*)(sbuf2 + nrqr); 842d5b485abSSatish Balay req_source = req_size + nrqr; 843d5b485abSSatish Balay 844d5b485abSSatish Balay { 845df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 846*b24ad042SBarry Smith PetscInt *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 847d5b485abSSatish Balay 848d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 849999d9058SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 850999d9058SBarry Smith req_size[idex] = 0; 851999d9058SBarry Smith rbuf1_i = rbuf1[idex]; 852d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 853*b24ad042SBarry Smith ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 854*b24ad042SBarry Smith ierr = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr); 855999d9058SBarry Smith sbuf2_i = sbuf2[idex]; 856d5b485abSSatish Balay for (j=start; j<end; j++) { 857d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 858d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 859d5b485abSSatish Balay sbuf2_i[j] = ncols; 860999d9058SBarry Smith req_size[idex] += ncols; 861d5b485abSSatish Balay } 862999d9058SBarry Smith req_source[idex] = r_status1[i].MPI_SOURCE; 863d5b485abSSatish Balay /* form the header */ 864999d9058SBarry Smith sbuf2_i[0] = req_size[idex]; 865d5b485abSSatish Balay for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; } 866*b24ad042SBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 867d5b485abSSatish Balay } 868d5b485abSSatish Balay } 869606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 870606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 871d5b485abSSatish Balay 872d5b485abSSatish Balay /* recv buffer sizes */ 873d5b485abSSatish Balay /* Receive messages*/ 874d5b485abSSatish Balay 875*b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr); 876b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr); 877b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr); 878b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr); 879b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr); 880d5b485abSSatish Balay 881d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 882999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 883*b24ad042SBarry Smith ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr); 884999d9058SBarry Smith ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr); 885*b24ad042SBarry Smith ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 886*b24ad042SBarry Smith ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 887d5b485abSSatish Balay } 888606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 889606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 890d5b485abSSatish Balay 891d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 892b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr); 893b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr); 894d5b485abSSatish Balay 895ca161407SBarry Smith ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr); 896ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr); 897606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 898606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 899606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 900606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 901d5b485abSSatish Balay 902d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 903*b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr); 904d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 905*b24ad042SBarry Smith ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr); 906d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 907d5b485abSSatish Balay 908b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr); 909d5b485abSSatish Balay { 910d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 911d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 912d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 913d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 914d5b485abSSatish Balay ct2 = 0; 915d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 916d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 917d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 918e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 919d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 920d5b485abSSatish Balay ncols = nzA + nzB; 921d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 922d5b485abSSatish Balay 923d5b485abSSatish Balay /* load the column indices for this row into cols*/ 924d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 925d5b485abSSatish Balay for (l=0; l<nzB; l++) { 926d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 927d5b485abSSatish Balay else break; 928d5b485abSSatish Balay } 929d5b485abSSatish Balay imark = l; 930d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 931d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 932d5b485abSSatish Balay ct2 += ncols; 933d5b485abSSatish Balay } 934d5b485abSSatish Balay } 935*b24ad042SBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 936d5b485abSSatish Balay } 937d5b485abSSatish Balay } 938b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr); 939b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr); 940d5b485abSSatish Balay 941d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 94282502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr); 943d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 94482502324SSatish Balay ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr); 945a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 946d5b485abSSatish Balay 947b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr); 948d5b485abSSatish Balay { 949d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 950d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 951d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 952d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 953d5b485abSSatish Balay ct2 = 0; 954d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 955d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 956d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 957e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 958d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 959d5b485abSSatish Balay ncols = nzA + nzB; 960d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 961a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 962d5b485abSSatish Balay 963d5b485abSSatish Balay /* load the column values for this row into vals*/ 9645b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 965d5b485abSSatish Balay for (l=0; l<nzB; l++) { 9663a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 9673eda8832SBarry Smith ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9683a40ed3dSBarry Smith } 969d5b485abSSatish Balay else break; 970d5b485abSSatish Balay } 971d5b485abSSatish Balay imark = l; 9723a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 9733eda8832SBarry Smith ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9743a40ed3dSBarry Smith } 9753a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 9763eda8832SBarry Smith ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9773a40ed3dSBarry Smith } 978d5b485abSSatish Balay ct2 += ncols; 979d5b485abSSatish Balay } 980d5b485abSSatish Balay } 9813eda8832SBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 982d5b485abSSatish Balay } 983d5b485abSSatish Balay } 984b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr); 985b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr); 986606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 987d5b485abSSatish Balay 988d5b485abSSatish Balay /* Form the matrix */ 989d5b485abSSatish Balay /* create col map */ 990d5b485abSSatish Balay { 991*b24ad042SBarry Smith PetscInt *icol_i; 992233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 993233c2ff6SSatish Balay /* Create row map*/ 994b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&colmaps);CHKERRQ(ierr); 995ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 996ff0794f7SSatish Balay ierr = PetscTableCreate(ncol[i]+1,&colmaps[i]);CHKERRQ(ierr); 997233c2ff6SSatish Balay } 998233c2ff6SSatish Balay #else 999*b24ad042SBarry Smith len = (1+ismax)*sizeof(PetscInt*)+ ismax*c->Nbs*sizeof(PetscInt); 100082502324SSatish Balay ierr = PetscMalloc(len,&cmap);CHKERRQ(ierr); 1001*b24ad042SBarry Smith cmap[0] = (PetscInt*)(cmap + ismax); 1002*b24ad042SBarry Smith ierr = PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(PetscInt));CHKERRQ(ierr); 1003a2feddc5SSatish Balay for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; } 1004233c2ff6SSatish Balay #endif 1005d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1006d5b485abSSatish Balay jmax = ncol[i]; 1007d5b485abSSatish Balay icol_i = icol[i]; 1008233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1009233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1010233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 1011001aedefSBarry Smith ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr); 1012233c2ff6SSatish Balay } 1013233c2ff6SSatish Balay #else 1014d5b485abSSatish Balay cmap_i = cmap[i]; 1015d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1016d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1017d5b485abSSatish Balay } 1018233c2ff6SSatish Balay #endif 1019d5b485abSSatish Balay } 1020d5b485abSSatish Balay } 1021d5b485abSSatish Balay 1022d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 1023d5b485abSSatish Balay for (i=0,j=0; i<ismax; i++) { j += nrow[i]; } 1024*b24ad042SBarry Smith len = (1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt); 102582502324SSatish Balay ierr = PetscMalloc(len,&lens);CHKERRQ(ierr); 1026*b24ad042SBarry Smith lens[0] = (PetscInt*)(lens + ismax); 1027*b24ad042SBarry Smith ierr = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr); 1028d5b485abSSatish Balay for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; } 1029d5b485abSSatish Balay 1030d5b485abSSatish Balay /* Update lens from local data */ 1031d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1032d5b485abSSatish Balay jmax = nrow[i]; 1033233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1034233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1035233c2ff6SSatish Balay #else 1036d5b485abSSatish Balay cmap_i = cmap[i]; 1037233c2ff6SSatish Balay #endif 1038d5b485abSSatish Balay irow_i = irow[i]; 1039d5b485abSSatish Balay lens_i = lens[i]; 1040d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1041d5b485abSSatish Balay row = irow_i[j]; 1042233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1043233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 1044233c2ff6SSatish Balay #else 1045d5b485abSSatish Balay proc = rtable[row]; 1046233c2ff6SSatish Balay #endif 1047d5b485abSSatish Balay if (proc == rank) { 104836f4e84dSSatish Balay /* Get indices from matA and then from matB */ 104936f4e84dSSatish Balay row = row - rstart; 105036f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 105136f4e84dSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1052233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1053233c2ff6SSatish Balay for (k=0; k<nzA; k++) { 1054233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1055233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1056233c2ff6SSatish Balay } 1057233c2ff6SSatish Balay for (k=0; k<nzB; k++) { 1058233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1059233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1060233c2ff6SSatish Balay } 1061233c2ff6SSatish Balay #else 1062ca161407SBarry Smith for (k=0; k<nzA; k++) { 106336f4e84dSSatish Balay if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; } 1064ca161407SBarry Smith } 1065ca161407SBarry Smith for (k=0; k<nzB; k++) { 106636f4e84dSSatish Balay if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; } 1067d5b485abSSatish Balay } 1068233c2ff6SSatish Balay #endif 1069d5b485abSSatish Balay } 1070a2feddc5SSatish Balay } 1071ca161407SBarry Smith } 1072233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1073d5b485abSSatish Balay /* Create row map*/ 1074b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rowmaps);CHKERRQ(ierr); 1075ff0794f7SSatish Balay for (i=0; i<ismax; i++){ 1076ff0794f7SSatish Balay ierr = PetscTableCreate(nrow[i]+1,&rowmaps[i]);CHKERRQ(ierr); 1077233c2ff6SSatish Balay } 1078233c2ff6SSatish Balay #else 1079233c2ff6SSatish Balay /* Create row map*/ 1080*b24ad042SBarry Smith len = (1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt); 108182502324SSatish Balay ierr = PetscMalloc(len,&rmap);CHKERRQ(ierr); 1082*b24ad042SBarry Smith rmap[0] = (PetscInt*)(rmap + ismax); 1083*b24ad042SBarry Smith ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr); 1084233c2ff6SSatish Balay for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;} 1085233c2ff6SSatish Balay #endif 1086d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1087d5b485abSSatish Balay irow_i = irow[i]; 1088d5b485abSSatish Balay jmax = nrow[i]; 1089233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1090233c2ff6SSatish Balay lrow1_grow1 = rowmaps[i]; 1091233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 1092233c2ff6SSatish Balay ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr); 1093233c2ff6SSatish Balay } 1094233c2ff6SSatish Balay #else 1095233c2ff6SSatish Balay rmap_i = rmap[i]; 1096d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1097d5b485abSSatish Balay rmap_i[irow_i[j]] = j; 1098d5b485abSSatish Balay } 1099233c2ff6SSatish Balay #endif 1100d5b485abSSatish Balay } 1101d5b485abSSatish Balay 1102d5b485abSSatish Balay /* Update lens from offproc data */ 1103d5b485abSSatish Balay { 1104*b24ad042SBarry Smith PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1105*b24ad042SBarry Smith PetscMPIInt ii; 1106d5b485abSSatish Balay 1107d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1108*b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr); 1109*b24ad042SBarry Smith idex = pa[ii]; 1110999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1111d5b485abSSatish Balay jmax = sbuf1_i[0]; 1112d5b485abSSatish Balay ct1 = 2*jmax+1; 1113d5b485abSSatish Balay ct2 = 0; 1114*b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1115*b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 1116d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1117d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1118d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1119d5b485abSSatish Balay lens_i = lens[is_no]; 1120233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1121233c2ff6SSatish Balay lcol1_gcol1 = colmaps[is_no]; 1122233c2ff6SSatish Balay lrow1_grow1 = rowmaps[is_no]; 1123233c2ff6SSatish Balay #else 1124d5b485abSSatish Balay cmap_i = cmap[is_no]; 1125d5b485abSSatish Balay rmap_i = rmap[is_no]; 1126233c2ff6SSatish Balay #endif 1127d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1128233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1129233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1130233c2ff6SSatish Balay row--; 1131e005ede5SBarry Smith if (row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); } 1132233c2ff6SSatish Balay #else 1133d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1134233c2ff6SSatish Balay #endif 1135d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1136d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1137233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1138233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 1139233c2ff6SSatish Balay if (tt) { 1140233c2ff6SSatish Balay lens_i[row]++; 1141233c2ff6SSatish Balay } 1142233c2ff6SSatish Balay #else 1143d5b485abSSatish Balay if (cmap_i[rbuf3_i[ct2]]) { 1144d5b485abSSatish Balay lens_i[row]++; 1145d5b485abSSatish Balay } 1146233c2ff6SSatish Balay #endif 1147d5b485abSSatish Balay } 1148d5b485abSSatish Balay } 1149d5b485abSSatish Balay } 1150d5b485abSSatish Balay } 1151d5b485abSSatish Balay } 1152606d414cSSatish Balay ierr = PetscFree(r_status3);CHKERRQ(ierr); 1153606d414cSSatish Balay ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1154ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr); 1155606d414cSSatish Balay ierr = PetscFree(s_status3);CHKERRQ(ierr); 1156606d414cSSatish Balay ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1157d5b485abSSatish Balay 1158d5b485abSSatish Balay /* Create the submatrices */ 1159d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 1160d5b485abSSatish Balay /* 1161d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1162d5b485abSSatish Balay */ 1163d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1164df36731bSBarry Smith mat = (Mat_SeqBAIJ *)(submats[i]->data); 116536f4e84dSSatish Balay if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || mat->bs != bs)) { 116629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 1167d5b485abSSatish Balay } 1168*b24ad042SBarry Smith ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 11690f5bd95cSBarry Smith if (flag == PETSC_FALSE) { 117029bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 1171d5b485abSSatish Balay } 1172d5b485abSSatish Balay /* Initial matrix as if empty */ 1173*b24ad042SBarry Smith ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr); 1174ce60de66SLois Curfman McInnes submats[i]->factor = C->factor; 1175d5b485abSSatish Balay } 1176ca161407SBarry Smith } else { 1177d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1178e7ef3d9dSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,nrow[i]*bs,ncol[i]*bs,nrow[i]*bs,ncol[i]*bs,submats+i);CHKERRQ(ierr); 1179e7ef3d9dSBarry Smith ierr = MatSetType(submats[i],A->type_name);CHKERRQ(ierr); 1180e7ef3d9dSBarry Smith ierr = MatSeqBAIJSetPreallocation(submats[i],a->bs,0,lens[i]);CHKERRQ(ierr); 1181*b24ad042SBarry Smith #if !defined(PETSC_USE_64BIT_INT) 1182e7ef3d9dSBarry Smith ierr = MatSeqSBAIJSetPreallocation(submats[i],a->bs,0,lens[i]);CHKERRQ(ierr); 1183*b24ad042SBarry Smith #endif 1184d5b485abSSatish Balay } 1185d5b485abSSatish Balay } 1186d5b485abSSatish Balay 1187d5b485abSSatish Balay /* Assemble the matrices */ 1188d5b485abSSatish Balay /* First assemble the local rows */ 1189d5b485abSSatish Balay { 1190*b24ad042SBarry Smith PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 11913eda8832SBarry Smith MatScalar *imat_a; 1192d5b485abSSatish Balay 1193d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1194df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[i]->data; 1195d5b485abSSatish Balay imat_ilen = mat->ilen; 1196d5b485abSSatish Balay imat_j = mat->j; 1197d5b485abSSatish Balay imat_i = mat->i; 1198d5b485abSSatish Balay imat_a = mat->a; 1199233c2ff6SSatish Balay 1200233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1201233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1202233c2ff6SSatish Balay lrow1_grow1 = rowmaps[i]; 1203233c2ff6SSatish Balay #else 1204d5b485abSSatish Balay cmap_i = cmap[i]; 1205d5b485abSSatish Balay rmap_i = rmap[i]; 1206233c2ff6SSatish Balay #endif 1207d5b485abSSatish Balay irow_i = irow[i]; 1208d5b485abSSatish Balay jmax = nrow[i]; 1209d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1210d5b485abSSatish Balay row = irow_i[j]; 1211233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1212233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 1213233c2ff6SSatish Balay #else 1214d5b485abSSatish Balay proc = rtable[row]; 1215233c2ff6SSatish Balay #endif 1216d5b485abSSatish Balay if (proc == rank) { 121736f4e84dSSatish Balay row = row - rstart; 121836f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 121936f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 122036f4e84dSSatish Balay cworkA = a_j + a_i[row]; 122136f4e84dSSatish Balay cworkB = b_j + b_i[row]; 122236f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 122336f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 1224233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1225233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr); 1226233c2ff6SSatish Balay row--; 1227e005ede5SBarry Smith if (row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); } 1228233c2ff6SSatish Balay #else 122936f4e84dSSatish Balay row = rmap_i[row + rstart]; 1230233c2ff6SSatish Balay #endif 1231df36731bSBarry Smith mat_i = imat_i[row]; 123236f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1233d5b485abSSatish Balay mat_j = imat_j + mat_i; 123436f4e84dSSatish Balay ilen_row = imat_ilen[row]; 123536f4e84dSSatish Balay 123636f4e84dSSatish Balay /* load the column indices for this row into cols*/ 123736f4e84dSSatish Balay for (l=0; l<nzB; l++) { 123836f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 1239233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1240233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr); 1241233c2ff6SSatish Balay if (tcol) { 1242233c2ff6SSatish Balay #else 124336f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1244233c2ff6SSatish Balay #endif 1245df36731bSBarry Smith *mat_j++ = tcol - 1; 12463eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1247549d3d68SSatish Balay mat_a += bs2; 1248d5b485abSSatish Balay ilen_row++; 1249d5b485abSSatish Balay } 1250ca161407SBarry Smith } else break; 125136f4e84dSSatish Balay } 125236f4e84dSSatish Balay imark = l; 125336f4e84dSSatish Balay for (l=0; l<nzA; l++) { 1254233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1255233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1256233c2ff6SSatish Balay if (tcol) { 1257233c2ff6SSatish Balay #else 125836f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 1259233c2ff6SSatish Balay #endif 126036f4e84dSSatish Balay *mat_j++ = tcol - 1; 12613eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1262549d3d68SSatish Balay mat_a += bs2; 126336f4e84dSSatish Balay ilen_row++; 126436f4e84dSSatish Balay } 126536f4e84dSSatish Balay } 126636f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 1267233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1268233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1269233c2ff6SSatish Balay if (tcol) { 1270233c2ff6SSatish Balay #else 127136f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1272233c2ff6SSatish Balay #endif 127336f4e84dSSatish Balay *mat_j++ = tcol - 1; 12743eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1275549d3d68SSatish Balay mat_a += bs2; 127636f4e84dSSatish Balay ilen_row++; 127736f4e84dSSatish Balay } 127836f4e84dSSatish Balay } 1279d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1280d5b485abSSatish Balay } 1281d5b485abSSatish Balay } 128236f4e84dSSatish Balay 1283d5b485abSSatish Balay } 1284d5b485abSSatish Balay } 1285d5b485abSSatish Balay 1286d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1287d5b485abSSatish Balay { 1288*b24ad042SBarry Smith PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1289*b24ad042SBarry Smith PetscInt *imat_j,*imat_i; 12903eda8832SBarry Smith MatScalar *imat_a,*rbuf4_i; 1291*b24ad042SBarry Smith PetscMPIInt ii; 1292d5b485abSSatish Balay 1293d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1294*b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr); 1295*b24ad042SBarry Smith idex = pa[ii]; 1296999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1297d5b485abSSatish Balay jmax = sbuf1_i[0]; 1298d5b485abSSatish Balay ct1 = 2*jmax + 1; 1299d5b485abSSatish Balay ct2 = 0; 1300*b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1301*b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 1302*b24ad042SBarry Smith rbuf4_i = rbuf4[ii]; 1303d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1304d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1305233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1306233c2ff6SSatish Balay lrow1_grow1 = rowmaps[is_no]; 1307233c2ff6SSatish Balay lcol1_gcol1 = colmaps[is_no]; 1308233c2ff6SSatish Balay #else 1309d5b485abSSatish Balay rmap_i = rmap[is_no]; 1310d5b485abSSatish Balay cmap_i = cmap[is_no]; 1311233c2ff6SSatish Balay #endif 1312df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[is_no]->data; 1313d5b485abSSatish Balay imat_ilen = mat->ilen; 1314d5b485abSSatish Balay imat_j = mat->j; 1315d5b485abSSatish Balay imat_i = mat->i; 1316d5b485abSSatish Balay imat_a = mat->a; 1317d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1318d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1319d5b485abSSatish Balay row = sbuf1_i[ct1]; 1320233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1321233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr); 1322233c2ff6SSatish Balay row--; 1323e005ede5SBarry Smith if(row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); } 1324233c2ff6SSatish Balay #else 1325d5b485abSSatish Balay row = rmap_i[row]; 1326233c2ff6SSatish Balay #endif 1327d5b485abSSatish Balay ilen = imat_ilen[row]; 1328df36731bSBarry Smith mat_i = imat_i[row]; 132936f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1330d5b485abSSatish Balay mat_j = imat_j + mat_i; 1331d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1332d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1333233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1334233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1335233c2ff6SSatish Balay if (tcol) { 1336233c2ff6SSatish Balay #else 1337d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1338233c2ff6SSatish Balay #endif 1339df36731bSBarry Smith *mat_j++ = tcol - 1; 134036f4e84dSSatish Balay /* *mat_a++= rbuf4_i[ct2]; */ 13413eda8832SBarry Smith ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1342549d3d68SSatish Balay mat_a += bs2; 1343d5b485abSSatish Balay ilen++; 1344d5b485abSSatish Balay } 1345d5b485abSSatish Balay } 1346d5b485abSSatish Balay imat_ilen[row] = ilen; 1347d5b485abSSatish Balay } 1348d5b485abSSatish Balay } 1349d5b485abSSatish Balay } 1350d5b485abSSatish Balay } 1351606d414cSSatish Balay ierr = PetscFree(r_status4);CHKERRQ(ierr); 1352606d414cSSatish Balay ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1353ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr); 1354606d414cSSatish Balay ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1355606d414cSSatish Balay ierr = PetscFree(s_status4);CHKERRQ(ierr); 1356d5b485abSSatish Balay 1357d5b485abSSatish Balay /* Restore the indices */ 1358d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1359d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1360d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1361d5b485abSSatish Balay } 1362d5b485abSSatish Balay 1363d5b485abSSatish Balay /* Destroy allocated memory */ 1364606d414cSSatish Balay ierr = PetscFree(irow);CHKERRQ(ierr); 1365606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 1366606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 1367d5b485abSSatish Balay 1368606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 1369606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 1370d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1371606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1372d5b485abSSatish Balay } 1373d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1374606d414cSSatish Balay ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1375606d414cSSatish Balay ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 1376d5b485abSSatish Balay } 1377d5b485abSSatish Balay 1378606d414cSSatish Balay ierr = PetscFree(sbuf2);CHKERRQ(ierr); 1379606d414cSSatish Balay ierr = PetscFree(rbuf3);CHKERRQ(ierr); 1380606d414cSSatish Balay ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1381606d414cSSatish Balay ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1382606d414cSSatish Balay ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1383606d414cSSatish Balay ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1384606d414cSSatish Balay ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1385d5b485abSSatish Balay 1386233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1387ff0794f7SSatish Balay for (i=0; i<ismax; i++){ 1388233c2ff6SSatish Balay ierr = PetscTableDelete(rowmaps[i]);CHKERRQ(ierr); 1389233c2ff6SSatish Balay ierr = PetscTableDelete(colmaps[i]);CHKERRQ(ierr); 1390233c2ff6SSatish Balay } 1391233c2ff6SSatish Balay ierr = PetscFree(colmaps);CHKERRQ(ierr); 1392233c2ff6SSatish Balay ierr = PetscFree(rowmaps);CHKERRQ(ierr); 1393233c2ff6SSatish Balay #else 1394606d414cSSatish Balay ierr = PetscFree(rmap);CHKERRQ(ierr); 1395233c2ff6SSatish Balay ierr = PetscFree(cmap);CHKERRQ(ierr); 1396233c2ff6SSatish Balay #endif 1397606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1398d5b485abSSatish Balay 1399d5b485abSSatish Balay for (i=0; i<ismax; i++) { 140036f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 140136f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1402d5b485abSSatish Balay } 140334e6ae18SSatish Balay 14043a40ed3dSBarry Smith PetscFunctionReturn(0); 1405d5b485abSSatish Balay } 1406ca161407SBarry Smith 1407