1*be1d678aSKris Buschelman #define PETSCMAT_DLL 2*be1d678aSKris Buschelman 3d5b485abSSatish Balay /* 4d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 5d5b485abSSatish Balay and to find submatrices that were shared across processors. 6d5b485abSSatish Balay */ 770f55243SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h" 80a835dfdSSatish Balay #include "petscbt.h" 9d5b485abSSatish Balay 10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat,PetscInt,IS *); 11b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char **,PetscInt*,PetscInt**); 12b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt **,PetscInt**,PetscInt*); 13b24ad042SBarry Smith EXTERN PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 14b24ad042SBarry Smith EXTERN PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 15d5b485abSSatish Balay 164a2ae208SSatish Balay #undef __FUNCT__ 174a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ" 18b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov) 19d5b485abSSatish Balay { 206849ba73SBarry Smith PetscErrorCode ierr; 21521d7252SBarry Smith PetscInt i,N=C->N, bs=C->bs; 2236f4e84dSSatish Balay IS *is_new; 2336f4e84dSSatish Balay 243a40ed3dSBarry Smith PetscFunctionBegin; 2582502324SSatish Balay ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr); 26df36731bSBarry Smith /* Convert the indices into block format */ 2727f478b1SHong Zhang ierr = ISCompressIndicesGeneral(N,bs,imax,is,is_new);CHKERRQ(ierr); 2829bbc08cSBarry Smith if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");} 29d5b485abSSatish Balay for (i=0; i<ov; ++i) { 3036f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 31d5b485abSSatish Balay } 32ca161407SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);} 3327f478b1SHong Zhang ierr = ISExpandIndicesGeneral(N,bs,imax,is_new,is);CHKERRQ(ierr); 34ca161407SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);} 35606d414cSSatish Balay ierr = PetscFree(is_new);CHKERRQ(ierr); 363a40ed3dSBarry Smith PetscFunctionReturn(0); 37d5b485abSSatish Balay } 38d5b485abSSatish Balay 39d5b485abSSatish Balay /* 40d5b485abSSatish Balay Sample message format: 41d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 420e9b0e7eSHong Zhang to index sets is[1], is[5] 43d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 44d5b485abSSatish Balay ----------- 45d5b485abSSatish Balay mesg [1] = 1 => is[1] 46d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 47d5b485abSSatish Balay ----------- 48d5b485abSSatish Balay mesg [5] = 5 => is[5] 49d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 50d5b485abSSatish Balay ----------- 51d5b485abSSatish Balay mesg [7] 520e9b0e7eSHong Zhang mesg [n] data(is[1]) 53d5b485abSSatish Balay ----------- 54d5b485abSSatish Balay mesg[n+1] 55d5b485abSSatish Balay mesg[m] data(is[5]) 56d5b485abSSatish Balay ----------- 57d5b485abSSatish Balay 58d5b485abSSatish Balay Notes: 59d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 60d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 61d5b485abSSatish Balay */ 624a2ae208SSatish Balay #undef __FUNCT__ 634a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once" 64b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[]) 65d5b485abSSatish Balay { 66df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 67b24ad042SBarry Smith PetscInt **idx,*n,*w3,*w4,*rtable,**data,len,*idx_i; 686849ba73SBarry Smith PetscErrorCode ierr; 69b24ad042SBarry Smith PetscMPIInt size,rank,tag1,tag2,*w2,*w1,nrqr; 70b24ad042SBarry Smith PetscInt Mbs,i,j,k,**rbuf,row,proc,nrqs,msz,**outdat,**ptr; 71b24ad042SBarry Smith PetscInt *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2; 72b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1,*onodes2,*olengths2; 73f1af5d2fSBarry Smith PetscBT *table; 74d5b485abSSatish Balay MPI_Comm comm; 75d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 76d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 77d5b485abSSatish Balay 783a40ed3dSBarry Smith PetscFunctionBegin; 79d5b485abSSatish Balay comm = 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 87b24ad042SBarry Smith len = (imax+1)*sizeof(PetscInt*)+ (imax + Mbs)*sizeof(PetscInt); 8882502324SSatish Balay ierr = PetscMalloc(len,&idx);CHKERRQ(ierr); 89b24ad042SBarry Smith n = (PetscInt*)(idx + imax); 90d5b485abSSatish Balay rtable = n + imax; 91d5b485abSSatish Balay 92d5b485abSSatish Balay for (i=0; i<imax; i++) { 93d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 94b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 95d5b485abSSatish Balay } 96d5b485abSSatish Balay 97d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 98d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 99d5b485abSSatish Balay len = c->rowners[i+1]; 100d5b485abSSatish Balay for (; j<len; j++) { 101d5b485abSSatish Balay rtable[j] = i; 102d5b485abSSatish Balay } 103d5b485abSSatish Balay } 104d5b485abSSatish Balay 105d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 106d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 107b24ad042SBarry Smith ierr = PetscMalloc(size*2*sizeof(PetscInt)+size*2*sizeof(PetscMPIInt),&w1);CHKERRQ(ierr);/* mesg size */ 108d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 109b24ad042SBarry Smith w3 = (PetscInt*) (w2 + size); /* no of IS that needs to be sent to proc i */ 110d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 111b24ad042SBarry Smith ierr = PetscMemzero(w1,size*sizeof(PetscInt)+2*size*sizeof(PetscMPIInt));CHKERRQ(ierr); /* initialise work vector*/ 112d5b485abSSatish Balay for (i=0; i<imax; i++) { 113b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 114d5b485abSSatish Balay idx_i = idx[i]; 115d5b485abSSatish Balay len = n[i]; 116d5b485abSSatish Balay for (j=0; j<len; j++) { 117d5b485abSSatish Balay row = idx_i[j]; 1186b41c64dSBarry Smith if (row < 0) { 119e005ede5SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 1206b41c64dSBarry Smith } 121d5b485abSSatish Balay proc = rtable[row]; 122d5b485abSSatish Balay w4[proc]++; 123d5b485abSSatish Balay } 124d5b485abSSatish Balay for (j=0; j<size; j++){ 125d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 126d5b485abSSatish Balay } 127d5b485abSSatish Balay } 128d5b485abSSatish Balay 129d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 130d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1310e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 132d5b485abSSatish Balay w3[rank] = 0; 133d5b485abSSatish Balay for (i=0; i<size; i++) { 134d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 135d5b485abSSatish Balay } 136d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 137b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); 138d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 139d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 140d5b485abSSatish Balay } 141d5b485abSSatish Balay 142d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 143d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 144d5b485abSSatish Balay j = pa[i]; 145d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 146d5b485abSSatish Balay msz += w1[j]; 147d5b485abSSatish Balay } 148d5b485abSSatish Balay 149c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 150a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 151c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 152d5b485abSSatish Balay 153c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 154c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 155d5b485abSSatish Balay 156d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 157b24ad042SBarry Smith len = 2*size*sizeof(PetscInt*) + (size+msz)*sizeof(PetscInt); 15882502324SSatish Balay ierr = PetscMalloc(len,&outdat);CHKERRQ(ierr); 159d5b485abSSatish Balay ptr = outdat + size; /* Pointers to the data in outgoing buffers */ 160b24ad042SBarry Smith ierr = PetscMemzero(outdat,2*size*sizeof(PetscInt*));CHKERRQ(ierr); 161b24ad042SBarry Smith tmp = (PetscInt*)(outdat + 2*size); 162d5b485abSSatish Balay ctr = tmp + msz; 163d5b485abSSatish Balay 164d5b485abSSatish Balay { 165b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 166d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 167d5b485abSSatish Balay j = pa[i]; 168d5b485abSSatish Balay iptr += ict; 169d5b485abSSatish Balay outdat[j] = iptr; 170d5b485abSSatish Balay ict = w1[j]; 171d5b485abSSatish Balay } 172d5b485abSSatish Balay } 173d5b485abSSatish Balay 174d5b485abSSatish Balay /* Form the outgoing messages */ 175d5b485abSSatish Balay /*plug in the headers*/ 176d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 177d5b485abSSatish Balay j = pa[i]; 178d5b485abSSatish Balay outdat[j][0] = 0; 179b24ad042SBarry Smith ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 180d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 181d5b485abSSatish Balay } 182d5b485abSSatish Balay 183d5b485abSSatish Balay /* Memory for doing local proc's work*/ 184d5b485abSSatish Balay { 185b24ad042SBarry Smith PetscInt *d_p; 186d5b485abSSatish Balay char *t_p; 187d5b485abSSatish Balay 188b24ad042SBarry Smith len = (imax)*(sizeof(PetscBT) + sizeof(PetscInt*)+ sizeof(PetscInt)) + 189b24ad042SBarry Smith (Mbs)*imax*sizeof(PetscInt) + (Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char) + 1; 19082502324SSatish Balay ierr = PetscMalloc(len,&table);CHKERRQ(ierr); 191549d3d68SSatish Balay ierr = PetscMemzero(table,len);CHKERRQ(ierr); 192b24ad042SBarry Smith data = (PetscInt **)(table + imax); 193b24ad042SBarry Smith isz = (PetscInt *)(data + imax); 194b24ad042SBarry Smith d_p = (PetscInt *)(isz + imax); 195bbd702dbSSatish Balay t_p = (char *)(d_p + Mbs*imax); 196bbd702dbSSatish Balay for (i=0; i<imax; i++) { 197b6410449SSatish Balay table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 198bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 199d5b485abSSatish Balay } 200d5b485abSSatish Balay } 201d5b485abSSatish Balay 202d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 203d5b485abSSatish Balay { 204b24ad042SBarry Smith PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 205f1af5d2fSBarry Smith PetscBT table_i; 206d5b485abSSatish Balay 207d5b485abSSatish Balay for (i=0; i<imax; i++) { 208b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 209d5b485abSSatish Balay n_i = n[i]; 210d5b485abSSatish Balay table_i = table[i]; 211d5b485abSSatish Balay idx_i = idx[i]; 212d5b485abSSatish Balay data_i = data[i]; 213d5b485abSSatish Balay isz_i = isz[i]; 214d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 215d5b485abSSatish Balay row = idx_i[j]; 216d5b485abSSatish Balay proc = rtable[row]; 217d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 218d5b485abSSatish Balay ctr[proc]++; 219d5b485abSSatish Balay *ptr[proc] = row; 220d5b485abSSatish Balay ptr[proc]++; 221d5b485abSSatish Balay } 222d5b485abSSatish Balay else { /* Update the local table */ 223f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 224d5b485abSSatish Balay } 225d5b485abSSatish Balay } 226d5b485abSSatish Balay /* Update the headers for the current IS */ 227d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 228d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 229d5b485abSSatish Balay outdat_j = outdat[j]; 230d5b485abSSatish Balay k = ++outdat_j[0]; 231d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 232d5b485abSSatish Balay outdat_j[2*k-1] = i; 233d5b485abSSatish Balay } 234d5b485abSSatish Balay } 235d5b485abSSatish Balay isz[i] = isz_i; 236d5b485abSSatish Balay } 237d5b485abSSatish Balay } 238d5b485abSSatish Balay 239d5b485abSSatish Balay /* Now post the sends */ 240b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 241d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 242d5b485abSSatish Balay j = pa[i]; 243b24ad042SBarry Smith ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 244d5b485abSSatish Balay } 245d5b485abSSatish Balay 246d5b485abSSatish Balay /* No longer need the original indices*/ 247d5b485abSSatish Balay for (i=0; i<imax; ++i) { 248d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 249d5b485abSSatish Balay } 250606d414cSSatish Balay ierr = PetscFree(idx);CHKERRQ(ierr); 251d5b485abSSatish Balay 252d5b485abSSatish Balay for (i=0; i<imax; ++i) { 253d5b485abSSatish Balay ierr = ISDestroy(is[i]);CHKERRQ(ierr); 254d5b485abSSatish Balay } 255d5b485abSSatish Balay 256d5b485abSSatish Balay /* Do Local work*/ 257df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 258d5b485abSSatish Balay 259d5b485abSSatish Balay /* Receive messages*/ 260b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr); 261c7dd2924SSatish Balay ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr); 262d5b485abSSatish Balay 263b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 264ca161407SBarry Smith ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr); 265d5b485abSSatish Balay 266d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 267606d414cSSatish Balay ierr = PetscFree(outdat);CHKERRQ(ierr); 268606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 269d5b485abSSatish Balay 270b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&xdata);CHKERRQ(ierr); 271b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt),&isz1);CHKERRQ(ierr); 272df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 273606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 274d5b485abSSatish Balay 275d5b485abSSatish Balay /* Send the data back*/ 276d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 277d5b485abSSatish Balay { 278b24ad042SBarry Smith PetscMPIInt *rw1; 279d5b485abSSatish Balay 280b24ad042SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&rw1);CHKERRQ(ierr); 281b24ad042SBarry Smith ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr); 282c7dd2924SSatish Balay 283d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 284d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 285e005ede5SBarry Smith if (proc != onodes1[i]) SETERRQ(PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 286d5b485abSSatish Balay rw1[proc] = isz1[i]; 287d5b485abSSatish Balay } 288d5b485abSSatish Balay 289c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 290c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 291d5b485abSSatish Balay 292c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 293c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 294c7dd2924SSatish Balay PetscFree(rw1); 295c7dd2924SSatish Balay } 296c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 297c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 298d5b485abSSatish Balay 299d5b485abSSatish Balay /* Now post the sends */ 300b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 301d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 302d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 303b24ad042SBarry Smith ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 304d5b485abSSatish Balay } 305d5b485abSSatish Balay 306d5b485abSSatish Balay /* receive work done on other processors*/ 307d5b485abSSatish Balay { 308b24ad042SBarry Smith PetscMPIInt idex; 309b24ad042SBarry Smith PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 310f1af5d2fSBarry Smith PetscBT table_i; 311d5b485abSSatish Balay MPI_Status *status2; 312d5b485abSSatish Balay 313169449a0SSatish Balay ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr); 314d5b485abSSatish Balay 315d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 316999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 317d5b485abSSatish Balay /* Process the message*/ 318999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 319d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 320999d9058SBarry Smith jmax = rbuf2[idex][0]; 321d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 322d5b485abSSatish Balay max = rbuf2_i[2*j]; 323d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 324d5b485abSSatish Balay isz_i = isz[is_no]; 325d5b485abSSatish Balay data_i = data[is_no]; 326d5b485abSSatish Balay table_i = table[is_no]; 327d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 328d5b485abSSatish Balay row = rbuf2_i[ct1]; 329f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 330d5b485abSSatish Balay } 331d5b485abSSatish Balay isz[is_no] = isz_i; 332d5b485abSSatish Balay } 333d5b485abSSatish Balay } 334ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr); 335606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 336d5b485abSSatish Balay } 337d5b485abSSatish Balay 338d5b485abSSatish Balay for (i=0; i<imax; ++i) { 339029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr); 340d5b485abSSatish Balay } 341d5b485abSSatish Balay 342c7dd2924SSatish Balay 343c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 344c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 345c7dd2924SSatish Balay 346606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 347606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 348606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 349606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 350606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 351606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 352606d414cSSatish Balay ierr = PetscFree(table);CHKERRQ(ierr); 353606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 354606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 355606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 356606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 357606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 3583a40ed3dSBarry Smith PetscFunctionReturn(0); 359d5b485abSSatish Balay } 360d5b485abSSatish Balay 3614a2ae208SSatish Balay #undef __FUNCT__ 3624a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local" 363d5b485abSSatish Balay /* 364df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 365d5b485abSSatish Balay the work on the local processor. 366d5b485abSSatish Balay 367d5b485abSSatish Balay Inputs: 368df36731bSBarry Smith C - MAT_MPIBAIJ; 369d5b485abSSatish Balay imax - total no of index sets processed at a time; 370df36731bSBarry Smith table - an array of char - size = Mbs bits. 371d5b485abSSatish Balay 372d5b485abSSatish Balay Output: 3730e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 374d5b485abSSatish Balay to each index set; 375d5b485abSSatish Balay data - pointer to the solutions 376d5b485abSSatish Balay */ 377b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 378d5b485abSSatish Balay { 379df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 380d5b485abSSatish Balay Mat A = c->A,B = c->B; 381df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 382b24ad042SBarry Smith PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 383b24ad042SBarry Smith PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 384f1af5d2fSBarry Smith PetscBT table_i; 385d5b485abSSatish Balay 3863a40ed3dSBarry Smith PetscFunctionBegin; 387d5b485abSSatish Balay rstart = c->rstart; 388d5b485abSSatish Balay cstart = c->cstart; 389d5b485abSSatish Balay ai = a->i; 390df36731bSBarry Smith aj = a->j; 391d5b485abSSatish Balay bi = b->i; 392df36731bSBarry Smith bj = b->j; 393d5b485abSSatish Balay garray = c->garray; 394d5b485abSSatish Balay 395d5b485abSSatish Balay 396d5b485abSSatish Balay for (i=0; i<imax; i++) { 397d5b485abSSatish Balay data_i = data[i]; 398d5b485abSSatish Balay table_i = table[i]; 399d5b485abSSatish Balay isz_i = isz[i]; 400d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 401d5b485abSSatish Balay row = data_i[j] - rstart; 402d5b485abSSatish Balay start = ai[row]; 403d5b485abSSatish Balay end = ai[row+1]; 404d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 405df36731bSBarry Smith val = aj[k] + cstart; 406f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 407d5b485abSSatish Balay } 408d5b485abSSatish Balay start = bi[row]; 409d5b485abSSatish Balay end = bi[row+1]; 410d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 411df36731bSBarry Smith val = garray[bj[k]]; 412f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 413d5b485abSSatish Balay } 414d5b485abSSatish Balay } 415d5b485abSSatish Balay isz[i] = isz_i; 416d5b485abSSatish Balay } 4173a40ed3dSBarry Smith PetscFunctionReturn(0); 418d5b485abSSatish Balay } 4194a2ae208SSatish Balay #undef __FUNCT__ 4204a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive" 421d5b485abSSatish Balay /* 422df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 423d5b485abSSatish Balay and return the output 424d5b485abSSatish Balay 425d5b485abSSatish Balay Input: 426d5b485abSSatish Balay C - the matrix 427d5b485abSSatish Balay nrqr - no of messages being processed. 428d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 429d5b485abSSatish Balay 430d5b485abSSatish Balay Output: 431d5b485abSSatish Balay xdata - array of messages to be sent back 432d5b485abSSatish Balay isz1 - size of each message 433d5b485abSSatish Balay 434d5b485abSSatish Balay For better efficiency perhaps we should malloc seperately each xdata[i], 435d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 4360e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of 437d5b485abSSatish Balay memory is used. 438d5b485abSSatish Balay 439d5b485abSSatish Balay */ 440b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 441d5b485abSSatish Balay { 442df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 443d5b485abSSatish Balay Mat A = c->A,B = c->B; 444df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 4456849ba73SBarry Smith PetscErrorCode ierr; 446b24ad042SBarry Smith PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 447b24ad042SBarry Smith PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 448b24ad042SBarry Smith PetscInt val,max1,max2,rank,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 449b24ad042SBarry Smith PetscInt *rbuf_i,kmax,rbuf_0; 450f1af5d2fSBarry Smith PetscBT xtable; 451d5b485abSSatish Balay 4523a40ed3dSBarry Smith PetscFunctionBegin; 453d5b485abSSatish Balay rank = c->rank; 454df36731bSBarry Smith Mbs = c->Mbs; 455d5b485abSSatish Balay rstart = c->rstart; 456d5b485abSSatish Balay cstart = c->cstart; 457d5b485abSSatish Balay ai = a->i; 458df36731bSBarry Smith aj = a->j; 459d5b485abSSatish Balay bi = b->i; 460df36731bSBarry Smith bj = b->j; 461d5b485abSSatish Balay garray = c->garray; 462d5b485abSSatish Balay 463d5b485abSSatish Balay 464d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 465d5b485abSSatish Balay rbuf_i = rbuf[i]; 466d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 467d5b485abSSatish Balay ct += rbuf_0; 468d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; } 469d5b485abSSatish Balay } 470d5b485abSSatish Balay 471701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 472701b8814SBarry Smith else max1 = 1; 473d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 474b24ad042SBarry Smith ierr = PetscMalloc(mem_estimate*sizeof(PetscInt),&xdata[0]);CHKERRQ(ierr); 475d5b485abSSatish Balay ++no_malloc; 4766831982aSBarry Smith ierr = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr); 477b24ad042SBarry Smith ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr); 478d5b485abSSatish Balay 479d5b485abSSatish Balay ct3 = 0; 480d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 481d5b485abSSatish Balay rbuf_i = rbuf[i]; 482d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 483d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 484d5b485abSSatish Balay ct2 = ct1; 485d5b485abSSatish Balay ct3 += ct1; 486d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 4876831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 488d5b485abSSatish Balay oct2 = ct2; 489d5b485abSSatish Balay kmax = rbuf_i[2*j]; 490d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 491d5b485abSSatish Balay row = rbuf_i[ct1]; 492f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 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++] = row; 503d5b485abSSatish Balay ct3++; 504d5b485abSSatish Balay } 505d5b485abSSatish Balay } 506d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 507d5b485abSSatish Balay row = xdata[i][k] - rstart; 508d5b485abSSatish Balay start = ai[row]; 509d5b485abSSatish Balay end = ai[row+1]; 510d5b485abSSatish Balay for (l=start; l<end; l++) { 511df36731bSBarry Smith val = aj[l] + cstart; 512f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 513d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 514b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 515b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 516b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 517606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 518d5b485abSSatish Balay xdata[0] = tmp; 519d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 520d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 521d5b485abSSatish Balay } 522d5b485abSSatish Balay xdata[i][ct2++] = val; 523d5b485abSSatish Balay ct3++; 524d5b485abSSatish Balay } 525d5b485abSSatish Balay } 526d5b485abSSatish Balay start = bi[row]; 527d5b485abSSatish Balay end = bi[row+1]; 528d5b485abSSatish Balay for (l=start; l<end; l++) { 529df36731bSBarry Smith val = garray[bj[l]]; 530f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 531d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 532b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 533b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 534b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 535606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 536d5b485abSSatish Balay xdata[0] = tmp; 537d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 538d5b485abSSatish Balay for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 539d5b485abSSatish Balay } 540d5b485abSSatish Balay xdata[i][ct2++] = val; 541d5b485abSSatish Balay ct3++; 542d5b485abSSatish Balay } 543d5b485abSSatish Balay } 544d5b485abSSatish Balay } 545d5b485abSSatish Balay /* Update the header*/ 546d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 547d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 548d5b485abSSatish Balay } 549d5b485abSSatish Balay xdata[i][0] = rbuf_0; 550d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 551d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 552d5b485abSSatish Balay } 5536831982aSBarry Smith ierr = PetscBTDestroy(xtable);CHKERRQ(ierr); 55463ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %D bytes, required %D, no of mallocs = %D\n",rank,mem_estimate,ct3,no_malloc));CHKERRQ(ierr); 5553a40ed3dSBarry Smith PetscFunctionReturn(0); 556d5b485abSSatish Balay } 557d5b485abSSatish Balay 558b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *); 559a2feddc5SSatish Balay 5604a2ae208SSatish Balay #undef __FUNCT__ 5614a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ" 562b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 563d5b485abSSatish Balay { 56436f4e84dSSatish Balay IS *isrow_new,*iscol_new; 565cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 5666849ba73SBarry Smith PetscErrorCode ierr; 567521d7252SBarry Smith PetscInt nmax,nstages_local,nstages,i,pos,max_no,N=C->N,bs=C->bs; 568a2feddc5SSatish Balay 5693a40ed3dSBarry Smith PetscFunctionBegin; 5703a40ed3dSBarry Smith /* The compression and expansion should be avoided. Does'nt point 571a2feddc5SSatish Balay out errors might change the indices hence buggey */ 572a2feddc5SSatish Balay 573b0a32e0cSBarry Smith ierr = PetscMalloc(2*(ismax+1)*sizeof(IS),&isrow_new);CHKERRQ(ierr); 57436f4e84dSSatish Balay iscol_new = isrow_new + ismax; 57527f478b1SHong Zhang ierr = ISCompressIndicesSorted(N,bs,ismax,isrow,isrow_new);CHKERRQ(ierr); 57627f478b1SHong Zhang ierr = ISCompressIndicesSorted(N,bs,ismax,iscol,iscol_new);CHKERRQ(ierr); 577cf4f527aSSatish Balay 578cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 579cf4f527aSSatish Balay if (scall != MAT_REUSE_MATRIX) { 58082502324SSatish Balay ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr); 581cf4f527aSSatish Balay } 582cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 583b24ad042SBarry Smith nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 584329f5518SBarry Smith if (!nmax) nmax = 1; 585cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 586cf4f527aSSatish Balay 587653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 588b24ad042SBarry Smith ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,C->comm);CHKERRQ(ierr); 589cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 590cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 591cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 592cf4f527aSSatish Balay else max_no = ismax-pos; 593cf4f527aSSatish Balay ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr); 594cf4f527aSSatish Balay pos += max_no; 595cf4f527aSSatish Balay } 59636f4e84dSSatish Balay 59736f4e84dSSatish Balay for (i=0; i<ismax; i++) { 598ca161407SBarry Smith ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr); 599ca161407SBarry Smith ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr); 60036f4e84dSSatish Balay } 601606d414cSSatish Balay ierr = PetscFree(isrow_new);CHKERRQ(ierr); 6023a40ed3dSBarry Smith PetscFunctionReturn(0); 603a2feddc5SSatish Balay } 604a2feddc5SSatish Balay 605233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 6064a2ae208SSatish Balay #undef __FUNCT__ 6074a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc" 608e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 609233c2ff6SSatish Balay { 610e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 611b24ad042SBarry Smith PetscMPIInt fproc = (PetscMPIInt) ((float)row * (float)size / (float)nGlobalNd + 0.5); 612233c2ff6SSatish Balay 613233c2ff6SSatish Balay PetscFunctionBegin; 61423ce1328SBarry Smith if (fproc > size) fproc = size; 615e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 616e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 617233c2ff6SSatish Balay else fproc++; 618233c2ff6SSatish Balay } 619e005ede5SBarry Smith *rank = fproc; 620233c2ff6SSatish Balay PetscFunctionReturn(0); 621233c2ff6SSatish Balay } 622233c2ff6SSatish Balay #endif 623233c2ff6SSatish Balay 624a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 625b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 6264a2ae208SSatish Balay #undef __FUNCT__ 6274a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local" 628b24ad042SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats) 629a2feddc5SSatish Balay { 630df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 631cf4f527aSSatish Balay Mat A = c->A; 632df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 633b24ad042SBarry Smith PetscInt **irow,**icol,*nrow,*ncol,*w3,*w4,start; 6346849ba73SBarry Smith PetscErrorCode ierr; 6359405f653SBarry Smith PetscMPIInt size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc; 6369405f653SBarry Smith PetscInt **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row; 637b24ad042SBarry Smith PetscInt nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol; 638b24ad042SBarry Smith PetscInt **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 639b24ad042SBarry Smith PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i; 640b24ad042SBarry Smith PetscInt len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 641521d7252SBarry Smith PetscInt bs=C->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 642b24ad042SBarry Smith PetscInt cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 643b24ad042SBarry Smith PetscInt *bmap = c->garray,ctmp,rstart=c->rstart; 644d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 645d5b485abSSatish Balay MPI_Request *r_waits4,*s_waits3,*s_waits4; 646d5b485abSSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 647d5b485abSSatish Balay MPI_Status *r_status3,*r_status4,*s_status4; 648d5b485abSSatish Balay MPI_Comm comm; 6493eda8832SBarry Smith MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB; 6503eda8832SBarry Smith MatScalar *a_a=a->a,*b_a=b->a; 6510f5bd95cSBarry Smith PetscTruth flag; 652b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1; 653c7dd2924SSatish Balay 65480d608c0SSatish Balay #if defined (PETSC_USE_CTABLE) 655b24ad042SBarry Smith PetscInt tt; 656233c2ff6SSatish Balay PetscTable *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1; 657233c2ff6SSatish Balay #else 658b24ad042SBarry Smith PetscInt **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 659233c2ff6SSatish Balay #endif 660d5b485abSSatish Balay 6613a40ed3dSBarry Smith PetscFunctionBegin; 662d5b485abSSatish Balay comm = C->comm; 66334e6ae18SSatish Balay tag0 = C->tag; 664d5b485abSSatish Balay size = c->size; 665d5b485abSSatish Balay rank = c->rank; 666d5b485abSSatish Balay 66734e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 66834e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 66934e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 67034e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 67134e6ae18SSatish Balay 672d5b485abSSatish Balay /* Check if the col indices are sorted */ 673d5b485abSSatish Balay for (i=0; i<ismax; i++) { 674888f2ed8SSatish Balay ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr); 67529bbc08cSBarry Smith if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted"); 676d5b485abSSatish Balay } 677d5b485abSSatish Balay 678b24ad042SBarry Smith len = (2*ismax+1)*(sizeof(PetscInt*)+ sizeof(PetscInt)); 679233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE) 680b24ad042SBarry Smith len += (Mbs+1)*sizeof(PetscInt); 681233c2ff6SSatish Balay #endif 68282502324SSatish Balay ierr = PetscMalloc(len,&irow);CHKERRQ(ierr); 683d5b485abSSatish Balay icol = irow + ismax; 684b24ad042SBarry Smith nrow = (PetscInt*)(icol + ismax); 685d5b485abSSatish Balay ncol = nrow + ismax; 686233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE) 687d5b485abSSatish Balay rtable = ncol + ismax; 688d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 689d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 690d5b485abSSatish Balay jmax = c->rowners[i+1]; 691d5b485abSSatish Balay for (; j<jmax; j++) { 692d5b485abSSatish Balay rtable[j] = i; 693d5b485abSSatish Balay } 694d5b485abSSatish Balay } 695233c2ff6SSatish Balay #endif 696233c2ff6SSatish Balay 697233c2ff6SSatish Balay for (i=0; i<ismax; i++) { 698233c2ff6SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 699233c2ff6SSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 700233c2ff6SSatish Balay ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 701233c2ff6SSatish Balay ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 702233c2ff6SSatish Balay } 703d5b485abSSatish Balay 704d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg,and buffer space 705d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 706b24ad042SBarry Smith ierr = PetscMalloc(size*2*sizeof(PetscInt) + size*2*sizeof(PetscMPIInt),&w1);CHKERRQ(ierr); /* mesg size */ 707d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 708b24ad042SBarry Smith w3 = (PetscInt*)(w2 + size); /* no of IS that needs to be sent to proc i */ 709d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 710b24ad042SBarry Smith ierr = PetscMemzero(w1,size*sizeof(PetscInt)+2*size*sizeof(PetscMPIInt));CHKERRQ(ierr); /* initialise work vector*/ 711d5b485abSSatish Balay for (i=0; i<ismax; i++) { 712b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 713d5b485abSSatish Balay jmax = nrow[i]; 714d5b485abSSatish Balay irow_i = irow[i]; 715d5b485abSSatish Balay for (j=0; j<jmax; j++) { 716d5b485abSSatish Balay row = irow_i[j]; 717233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 718233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 719233c2ff6SSatish Balay #else 720d5b485abSSatish Balay proc = rtable[row]; 721233c2ff6SSatish Balay #endif 722d5b485abSSatish Balay w4[proc]++; 723d5b485abSSatish Balay } 724d5b485abSSatish Balay for (j=0; j<size; j++) { 725d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 726d5b485abSSatish Balay } 727d5b485abSSatish Balay } 728d5b485abSSatish Balay 729d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 730e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 731d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 732d5b485abSSatish Balay w3[rank] = 0; 733d5b485abSSatish Balay for (i=0; i<size; i++) { 734d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 735d5b485abSSatish Balay } 736b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/ 737d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 738d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 739d5b485abSSatish Balay } 740d5b485abSSatish Balay 741d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 742d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 743d5b485abSSatish Balay j = pa[i]; 744d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 745d5b485abSSatish Balay msz += w1[j]; 746d5b485abSSatish Balay } 747d5b485abSSatish Balay 748c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 749a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 750c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 751d5b485abSSatish Balay 752c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 753c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 754c7dd2924SSatish Balay 755c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 756c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 757d5b485abSSatish Balay 758d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 759b24ad042SBarry Smith len = 2*size*sizeof(PetscInt*) + 2*msz*sizeof(PetscInt) + size*sizeof(PetscInt); 76082502324SSatish Balay ierr = PetscMalloc(len,&sbuf1);CHKERRQ(ierr); 761d5b485abSSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 762b24ad042SBarry Smith ierr = PetscMemzero(sbuf1,2*size*sizeof(PetscInt*));CHKERRQ(ierr); 763d5b485abSSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 764b24ad042SBarry Smith tmp = (PetscInt*)(ptr + size); 765d5b485abSSatish Balay ctr = tmp + 2*msz; 766d5b485abSSatish Balay 767d5b485abSSatish Balay { 768b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 769d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 770d5b485abSSatish Balay j = pa[i]; 771d5b485abSSatish Balay iptr += ict; 772d5b485abSSatish Balay sbuf1[j] = iptr; 773d5b485abSSatish Balay ict = w1[j]; 774d5b485abSSatish Balay } 775d5b485abSSatish Balay } 776d5b485abSSatish Balay 777d5b485abSSatish Balay /* Form the outgoing messages */ 778d5b485abSSatish Balay /* Initialise the header space */ 779d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 780d5b485abSSatish Balay j = pa[i]; 781d5b485abSSatish Balay sbuf1[j][0] = 0; 782b24ad042SBarry Smith ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 783d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 784d5b485abSSatish Balay } 785d5b485abSSatish Balay 786d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 787d5b485abSSatish Balay for (i=0; i<ismax; i++) { 788b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 789d5b485abSSatish Balay irow_i = irow[i]; 790d5b485abSSatish Balay jmax = nrow[i]; 791d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 792d5b485abSSatish Balay row = irow_i[j]; 793233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 794233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 795233c2ff6SSatish Balay #else 796d5b485abSSatish Balay proc = rtable[row]; 797233c2ff6SSatish Balay #endif 798d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 799d5b485abSSatish Balay ctr[proc]++; 800d5b485abSSatish Balay *ptr[proc] = row; 801d5b485abSSatish Balay ptr[proc]++; 802d5b485abSSatish Balay } 803d5b485abSSatish Balay } 804d5b485abSSatish Balay /* Update the headers for the current IS */ 805d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 80606ef35abSSatish Balay if ((ctr_j = ctr[j])) { 807d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 808d5b485abSSatish Balay k = ++sbuf1_j[0]; 809d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 810d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 811d5b485abSSatish Balay } 812d5b485abSSatish Balay } 813d5b485abSSatish Balay } 814d5b485abSSatish Balay 815d5b485abSSatish Balay /* Now post the sends */ 816b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 817d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 818d5b485abSSatish Balay j = pa[i]; 819b24ad042SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 820d5b485abSSatish Balay } 821d5b485abSSatish Balay 822d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 823b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr); 824b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr); 825d5b485abSSatish Balay rbuf2[0] = tmp + msz; 826d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 827d5b485abSSatish Balay j = pa[i]; 828d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 829d5b485abSSatish Balay } 830d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 831d5b485abSSatish Balay j = pa[i]; 832b24ad042SBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 833d5b485abSSatish Balay } 834d5b485abSSatish Balay 835d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 836d5b485abSSatish Balay 837d5b485abSSatish Balay /* Receive messages*/ 838b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 839b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr); 840b24ad042SBarry Smith len = 2*nrqr*sizeof(PetscInt) + (nrqr+1)*sizeof(PetscInt*); 84182502324SSatish Balay ierr = PetscMalloc(len,&sbuf2);CHKERRQ(ierr); 842b24ad042SBarry Smith req_size = (PetscInt*)(sbuf2 + nrqr); 843d5b485abSSatish Balay req_source = req_size + nrqr; 844d5b485abSSatish Balay 845d5b485abSSatish Balay { 846df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 847b24ad042SBarry Smith PetscInt *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 848d5b485abSSatish Balay 849d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 850999d9058SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 851999d9058SBarry Smith req_size[idex] = 0; 852999d9058SBarry Smith rbuf1_i = rbuf1[idex]; 853d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 854b24ad042SBarry Smith ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 855b24ad042SBarry Smith ierr = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr); 856999d9058SBarry Smith sbuf2_i = sbuf2[idex]; 857d5b485abSSatish Balay for (j=start; j<end; j++) { 858d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 859d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 860d5b485abSSatish Balay sbuf2_i[j] = ncols; 861999d9058SBarry Smith req_size[idex] += ncols; 862d5b485abSSatish Balay } 863999d9058SBarry Smith req_source[idex] = r_status1[i].MPI_SOURCE; 864d5b485abSSatish Balay /* form the header */ 865999d9058SBarry Smith sbuf2_i[0] = req_size[idex]; 866d5b485abSSatish Balay for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; } 867b24ad042SBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 868d5b485abSSatish Balay } 869d5b485abSSatish Balay } 870606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 871606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 872d5b485abSSatish Balay 873d5b485abSSatish Balay /* recv buffer sizes */ 874d5b485abSSatish Balay /* Receive messages*/ 875d5b485abSSatish Balay 876b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr); 877b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr); 878b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr); 879b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr); 880b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr); 881d5b485abSSatish Balay 882d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 883999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 884b24ad042SBarry Smith ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr); 885999d9058SBarry Smith ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr); 886b24ad042SBarry Smith ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 887b24ad042SBarry Smith ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 888d5b485abSSatish Balay } 889606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 890606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 891d5b485abSSatish Balay 892d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 893b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr); 894b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr); 895d5b485abSSatish Balay 896ca161407SBarry Smith ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr); 897ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr); 898606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 899606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 900606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 901606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 902d5b485abSSatish Balay 903d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 904b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr); 905d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 906b24ad042SBarry Smith ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr); 907d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 908d5b485abSSatish Balay 909b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr); 910d5b485abSSatish Balay { 911d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 912d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 913d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 914d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 915d5b485abSSatish Balay ct2 = 0; 916d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 917d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 918d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 919e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 920d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 921d5b485abSSatish Balay ncols = nzA + nzB; 922d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 923d5b485abSSatish Balay 924d5b485abSSatish Balay /* load the column indices for this row into cols*/ 925d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 926d5b485abSSatish Balay for (l=0; l<nzB; l++) { 927d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 928d5b485abSSatish Balay else break; 929d5b485abSSatish Balay } 930d5b485abSSatish Balay imark = l; 931d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 932d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 933d5b485abSSatish Balay ct2 += ncols; 934d5b485abSSatish Balay } 935d5b485abSSatish Balay } 936b24ad042SBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 937d5b485abSSatish Balay } 938d5b485abSSatish Balay } 939b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr); 940b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr); 941d5b485abSSatish Balay 942d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 94382502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr); 944d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 94582502324SSatish Balay ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr); 946a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 947d5b485abSSatish Balay 948b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr); 949d5b485abSSatish Balay { 950d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 951d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 952d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 953d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 954d5b485abSSatish Balay ct2 = 0; 955d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 956d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 957d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 958e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 959d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 960d5b485abSSatish Balay ncols = nzA + nzB; 961d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 962a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 963d5b485abSSatish Balay 964d5b485abSSatish Balay /* load the column values for this row into vals*/ 9655b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 966d5b485abSSatish Balay for (l=0; l<nzB; l++) { 9673a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 9683eda8832SBarry Smith ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9693a40ed3dSBarry Smith } 970d5b485abSSatish Balay else break; 971d5b485abSSatish Balay } 972d5b485abSSatish Balay imark = l; 9733a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 9743eda8832SBarry Smith ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9753a40ed3dSBarry Smith } 9763a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 9773eda8832SBarry Smith ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9783a40ed3dSBarry Smith } 979d5b485abSSatish Balay ct2 += ncols; 980d5b485abSSatish Balay } 981d5b485abSSatish Balay } 9823eda8832SBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 983d5b485abSSatish Balay } 984d5b485abSSatish Balay } 985b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr); 986b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr); 987606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 988d5b485abSSatish Balay 989d5b485abSSatish Balay /* Form the matrix */ 990d5b485abSSatish Balay /* create col map */ 991d5b485abSSatish Balay { 992b24ad042SBarry Smith PetscInt *icol_i; 993233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 994233c2ff6SSatish Balay /* Create row map*/ 995b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&colmaps);CHKERRQ(ierr); 996ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 997ff0794f7SSatish Balay ierr = PetscTableCreate(ncol[i]+1,&colmaps[i]);CHKERRQ(ierr); 998233c2ff6SSatish Balay } 999233c2ff6SSatish Balay #else 1000b24ad042SBarry Smith len = (1+ismax)*sizeof(PetscInt*)+ ismax*c->Nbs*sizeof(PetscInt); 100182502324SSatish Balay ierr = PetscMalloc(len,&cmap);CHKERRQ(ierr); 1002b24ad042SBarry Smith cmap[0] = (PetscInt*)(cmap + ismax); 1003b24ad042SBarry Smith ierr = PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(PetscInt));CHKERRQ(ierr); 1004a2feddc5SSatish Balay for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; } 1005233c2ff6SSatish Balay #endif 1006d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1007d5b485abSSatish Balay jmax = ncol[i]; 1008d5b485abSSatish Balay icol_i = icol[i]; 1009233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1010233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1011233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 1012001aedefSBarry Smith ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr); 1013233c2ff6SSatish Balay } 1014233c2ff6SSatish Balay #else 1015d5b485abSSatish Balay cmap_i = cmap[i]; 1016d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1017d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1018d5b485abSSatish Balay } 1019233c2ff6SSatish Balay #endif 1020d5b485abSSatish Balay } 1021d5b485abSSatish Balay } 1022d5b485abSSatish Balay 1023d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 1024d5b485abSSatish Balay for (i=0,j=0; i<ismax; i++) { j += nrow[i]; } 1025b24ad042SBarry Smith len = (1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt); 102682502324SSatish Balay ierr = PetscMalloc(len,&lens);CHKERRQ(ierr); 1027b24ad042SBarry Smith lens[0] = (PetscInt*)(lens + ismax); 1028b24ad042SBarry Smith ierr = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr); 1029d5b485abSSatish Balay for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; } 1030d5b485abSSatish Balay 1031d5b485abSSatish Balay /* Update lens from local data */ 1032d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1033d5b485abSSatish Balay jmax = nrow[i]; 1034233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1035233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1036233c2ff6SSatish Balay #else 1037d5b485abSSatish Balay cmap_i = cmap[i]; 1038233c2ff6SSatish Balay #endif 1039d5b485abSSatish Balay irow_i = irow[i]; 1040d5b485abSSatish Balay lens_i = lens[i]; 1041d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1042d5b485abSSatish Balay row = irow_i[j]; 1043233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1044233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 1045233c2ff6SSatish Balay #else 1046d5b485abSSatish Balay proc = rtable[row]; 1047233c2ff6SSatish Balay #endif 1048d5b485abSSatish Balay if (proc == rank) { 104936f4e84dSSatish Balay /* Get indices from matA and then from matB */ 105036f4e84dSSatish Balay row = row - rstart; 105136f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 105236f4e84dSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1053233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1054233c2ff6SSatish Balay for (k=0; k<nzA; k++) { 1055233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1056233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1057233c2ff6SSatish Balay } 1058233c2ff6SSatish Balay for (k=0; k<nzB; k++) { 1059233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1060233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1061233c2ff6SSatish Balay } 1062233c2ff6SSatish Balay #else 1063ca161407SBarry Smith for (k=0; k<nzA; k++) { 106436f4e84dSSatish Balay if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; } 1065ca161407SBarry Smith } 1066ca161407SBarry Smith for (k=0; k<nzB; k++) { 106736f4e84dSSatish Balay if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; } 1068d5b485abSSatish Balay } 1069233c2ff6SSatish Balay #endif 1070d5b485abSSatish Balay } 1071a2feddc5SSatish Balay } 1072ca161407SBarry Smith } 1073233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1074d5b485abSSatish Balay /* Create row map*/ 1075b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rowmaps);CHKERRQ(ierr); 1076ff0794f7SSatish Balay for (i=0; i<ismax; i++){ 1077ff0794f7SSatish Balay ierr = PetscTableCreate(nrow[i]+1,&rowmaps[i]);CHKERRQ(ierr); 1078233c2ff6SSatish Balay } 1079233c2ff6SSatish Balay #else 1080233c2ff6SSatish Balay /* Create row map*/ 1081b24ad042SBarry Smith len = (1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt); 108282502324SSatish Balay ierr = PetscMalloc(len,&rmap);CHKERRQ(ierr); 1083b24ad042SBarry Smith rmap[0] = (PetscInt*)(rmap + ismax); 1084b24ad042SBarry Smith ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr); 1085233c2ff6SSatish Balay for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;} 1086233c2ff6SSatish Balay #endif 1087d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1088d5b485abSSatish Balay irow_i = irow[i]; 1089d5b485abSSatish Balay jmax = nrow[i]; 1090233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1091233c2ff6SSatish Balay lrow1_grow1 = rowmaps[i]; 1092233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 1093233c2ff6SSatish Balay ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr); 1094233c2ff6SSatish Balay } 1095233c2ff6SSatish Balay #else 1096233c2ff6SSatish Balay rmap_i = rmap[i]; 1097d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1098d5b485abSSatish Balay rmap_i[irow_i[j]] = j; 1099d5b485abSSatish Balay } 1100233c2ff6SSatish Balay #endif 1101d5b485abSSatish Balay } 1102d5b485abSSatish Balay 1103d5b485abSSatish Balay /* Update lens from offproc data */ 1104d5b485abSSatish Balay { 1105b24ad042SBarry Smith PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1106b24ad042SBarry Smith PetscMPIInt ii; 1107d5b485abSSatish Balay 1108d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1109b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr); 1110b24ad042SBarry Smith idex = pa[ii]; 1111999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1112d5b485abSSatish Balay jmax = sbuf1_i[0]; 1113d5b485abSSatish Balay ct1 = 2*jmax+1; 1114d5b485abSSatish Balay ct2 = 0; 1115b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1116b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 1117d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1118d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1119d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1120d5b485abSSatish Balay lens_i = lens[is_no]; 1121233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1122233c2ff6SSatish Balay lcol1_gcol1 = colmaps[is_no]; 1123233c2ff6SSatish Balay lrow1_grow1 = rowmaps[is_no]; 1124233c2ff6SSatish Balay #else 1125d5b485abSSatish Balay cmap_i = cmap[is_no]; 1126d5b485abSSatish Balay rmap_i = rmap[is_no]; 1127233c2ff6SSatish Balay #endif 1128d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1129233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1130233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1131233c2ff6SSatish Balay row--; 1132e005ede5SBarry Smith if (row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); } 1133233c2ff6SSatish Balay #else 1134d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1135233c2ff6SSatish Balay #endif 1136d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1137d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1138233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1139233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 1140233c2ff6SSatish Balay if (tt) { 1141233c2ff6SSatish Balay lens_i[row]++; 1142233c2ff6SSatish Balay } 1143233c2ff6SSatish Balay #else 1144d5b485abSSatish Balay if (cmap_i[rbuf3_i[ct2]]) { 1145d5b485abSSatish Balay lens_i[row]++; 1146d5b485abSSatish Balay } 1147233c2ff6SSatish Balay #endif 1148d5b485abSSatish Balay } 1149d5b485abSSatish Balay } 1150d5b485abSSatish Balay } 1151d5b485abSSatish Balay } 1152d5b485abSSatish Balay } 1153606d414cSSatish Balay ierr = PetscFree(r_status3);CHKERRQ(ierr); 1154606d414cSSatish Balay ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1155ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr); 1156606d414cSSatish Balay ierr = PetscFree(s_status3);CHKERRQ(ierr); 1157606d414cSSatish Balay ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1158d5b485abSSatish Balay 1159d5b485abSSatish Balay /* Create the submatrices */ 1160d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 1161d5b485abSSatish Balay /* 1162d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1163d5b485abSSatish Balay */ 1164d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1165df36731bSBarry Smith mat = (Mat_SeqBAIJ *)(submats[i]->data); 1166521d7252SBarry Smith if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->bs != bs)) { 116729bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 1168d5b485abSSatish Balay } 1169b24ad042SBarry Smith ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 1170abc0a331SBarry Smith if (!flag) { 117129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 1172d5b485abSSatish Balay } 1173d5b485abSSatish Balay /* Initial matrix as if empty */ 1174b24ad042SBarry Smith ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr); 1175ce60de66SLois Curfman McInnes submats[i]->factor = C->factor; 1176d5b485abSSatish Balay } 1177ca161407SBarry Smith } else { 1178d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1179e7ef3d9dSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,nrow[i]*bs,ncol[i]*bs,nrow[i]*bs,ncol[i]*bs,submats+i);CHKERRQ(ierr); 1180e7ef3d9dSBarry Smith ierr = MatSetType(submats[i],A->type_name);CHKERRQ(ierr); 1181521d7252SBarry Smith ierr = MatSeqBAIJSetPreallocation(submats[i],C->bs,0,lens[i]);CHKERRQ(ierr); 1182b24ad042SBarry Smith #if !defined(PETSC_USE_64BIT_INT) 1183521d7252SBarry Smith ierr = MatSeqSBAIJSetPreallocation(submats[i],C->bs,0,lens[i]);CHKERRQ(ierr); 1184b24ad042SBarry Smith #endif 1185d5b485abSSatish Balay } 1186d5b485abSSatish Balay } 1187d5b485abSSatish Balay 1188d5b485abSSatish Balay /* Assemble the matrices */ 1189d5b485abSSatish Balay /* First assemble the local rows */ 1190d5b485abSSatish Balay { 1191b24ad042SBarry Smith PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 11923eda8832SBarry Smith MatScalar *imat_a; 1193d5b485abSSatish Balay 1194d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1195df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[i]->data; 1196d5b485abSSatish Balay imat_ilen = mat->ilen; 1197d5b485abSSatish Balay imat_j = mat->j; 1198d5b485abSSatish Balay imat_i = mat->i; 1199d5b485abSSatish Balay imat_a = mat->a; 1200233c2ff6SSatish Balay 1201233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1202233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1203233c2ff6SSatish Balay lrow1_grow1 = rowmaps[i]; 1204233c2ff6SSatish Balay #else 1205d5b485abSSatish Balay cmap_i = cmap[i]; 1206d5b485abSSatish Balay rmap_i = rmap[i]; 1207233c2ff6SSatish Balay #endif 1208d5b485abSSatish Balay irow_i = irow[i]; 1209d5b485abSSatish Balay jmax = nrow[i]; 1210d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1211d5b485abSSatish Balay row = irow_i[j]; 1212233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1213233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 1214233c2ff6SSatish Balay #else 1215d5b485abSSatish Balay proc = rtable[row]; 1216233c2ff6SSatish Balay #endif 1217d5b485abSSatish Balay if (proc == rank) { 121836f4e84dSSatish Balay row = row - rstart; 121936f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 122036f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 122136f4e84dSSatish Balay cworkA = a_j + a_i[row]; 122236f4e84dSSatish Balay cworkB = b_j + b_i[row]; 122336f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 122436f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 1225233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1226233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr); 1227233c2ff6SSatish Balay row--; 1228e005ede5SBarry Smith if (row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); } 1229233c2ff6SSatish Balay #else 123036f4e84dSSatish Balay row = rmap_i[row + rstart]; 1231233c2ff6SSatish Balay #endif 1232df36731bSBarry Smith mat_i = imat_i[row]; 123336f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1234d5b485abSSatish Balay mat_j = imat_j + mat_i; 123536f4e84dSSatish Balay ilen_row = imat_ilen[row]; 123636f4e84dSSatish Balay 123736f4e84dSSatish Balay /* load the column indices for this row into cols*/ 123836f4e84dSSatish Balay for (l=0; l<nzB; l++) { 123936f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 1240233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1241233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr); 1242233c2ff6SSatish Balay if (tcol) { 1243233c2ff6SSatish Balay #else 124436f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1245233c2ff6SSatish Balay #endif 1246df36731bSBarry Smith *mat_j++ = tcol - 1; 12473eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1248549d3d68SSatish Balay mat_a += bs2; 1249d5b485abSSatish Balay ilen_row++; 1250d5b485abSSatish Balay } 1251ca161407SBarry Smith } else break; 125236f4e84dSSatish Balay } 125336f4e84dSSatish Balay imark = l; 125436f4e84dSSatish Balay for (l=0; l<nzA; l++) { 1255233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1256233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1257233c2ff6SSatish Balay if (tcol) { 1258233c2ff6SSatish Balay #else 125936f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 1260233c2ff6SSatish Balay #endif 126136f4e84dSSatish Balay *mat_j++ = tcol - 1; 12623eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1263549d3d68SSatish Balay mat_a += bs2; 126436f4e84dSSatish Balay ilen_row++; 126536f4e84dSSatish Balay } 126636f4e84dSSatish Balay } 126736f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 1268233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1269233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1270233c2ff6SSatish Balay if (tcol) { 1271233c2ff6SSatish Balay #else 127236f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1273233c2ff6SSatish Balay #endif 127436f4e84dSSatish Balay *mat_j++ = tcol - 1; 12753eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1276549d3d68SSatish Balay mat_a += bs2; 127736f4e84dSSatish Balay ilen_row++; 127836f4e84dSSatish Balay } 127936f4e84dSSatish Balay } 1280d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1281d5b485abSSatish Balay } 1282d5b485abSSatish Balay } 128336f4e84dSSatish Balay 1284d5b485abSSatish Balay } 1285d5b485abSSatish Balay } 1286d5b485abSSatish Balay 1287d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1288d5b485abSSatish Balay { 1289b24ad042SBarry Smith PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1290b24ad042SBarry Smith PetscInt *imat_j,*imat_i; 12913eda8832SBarry Smith MatScalar *imat_a,*rbuf4_i; 1292b24ad042SBarry Smith PetscMPIInt ii; 1293d5b485abSSatish Balay 1294d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1295b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr); 1296b24ad042SBarry Smith idex = pa[ii]; 1297999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1298d5b485abSSatish Balay jmax = sbuf1_i[0]; 1299d5b485abSSatish Balay ct1 = 2*jmax + 1; 1300d5b485abSSatish Balay ct2 = 0; 1301b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1302b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 1303b24ad042SBarry Smith rbuf4_i = rbuf4[ii]; 1304d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1305d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1306233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1307233c2ff6SSatish Balay lrow1_grow1 = rowmaps[is_no]; 1308233c2ff6SSatish Balay lcol1_gcol1 = colmaps[is_no]; 1309233c2ff6SSatish Balay #else 1310d5b485abSSatish Balay rmap_i = rmap[is_no]; 1311d5b485abSSatish Balay cmap_i = cmap[is_no]; 1312233c2ff6SSatish Balay #endif 1313df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[is_no]->data; 1314d5b485abSSatish Balay imat_ilen = mat->ilen; 1315d5b485abSSatish Balay imat_j = mat->j; 1316d5b485abSSatish Balay imat_i = mat->i; 1317d5b485abSSatish Balay imat_a = mat->a; 1318d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1319d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1320d5b485abSSatish Balay row = sbuf1_i[ct1]; 1321233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1322233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr); 1323233c2ff6SSatish Balay row--; 1324e005ede5SBarry Smith if(row < 0) { SETERRQ(PETSC_ERR_PLIB,"row not found in table"); } 1325233c2ff6SSatish Balay #else 1326d5b485abSSatish Balay row = rmap_i[row]; 1327233c2ff6SSatish Balay #endif 1328d5b485abSSatish Balay ilen = imat_ilen[row]; 1329df36731bSBarry Smith mat_i = imat_i[row]; 133036f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1331d5b485abSSatish Balay mat_j = imat_j + mat_i; 1332d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1333d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1334233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1335233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1336233c2ff6SSatish Balay if (tcol) { 1337233c2ff6SSatish Balay #else 1338d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1339233c2ff6SSatish Balay #endif 1340df36731bSBarry Smith *mat_j++ = tcol - 1; 134136f4e84dSSatish Balay /* *mat_a++= rbuf4_i[ct2]; */ 13423eda8832SBarry Smith ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1343549d3d68SSatish Balay mat_a += bs2; 1344d5b485abSSatish Balay ilen++; 1345d5b485abSSatish Balay } 1346d5b485abSSatish Balay } 1347d5b485abSSatish Balay imat_ilen[row] = ilen; 1348d5b485abSSatish Balay } 1349d5b485abSSatish Balay } 1350d5b485abSSatish Balay } 1351d5b485abSSatish Balay } 1352606d414cSSatish Balay ierr = PetscFree(r_status4);CHKERRQ(ierr); 1353606d414cSSatish Balay ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1354ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr); 1355606d414cSSatish Balay ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1356606d414cSSatish Balay ierr = PetscFree(s_status4);CHKERRQ(ierr); 1357d5b485abSSatish Balay 1358d5b485abSSatish Balay /* Restore the indices */ 1359d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1360d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1361d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1362d5b485abSSatish Balay } 1363d5b485abSSatish Balay 1364d5b485abSSatish Balay /* Destroy allocated memory */ 1365606d414cSSatish Balay ierr = PetscFree(irow);CHKERRQ(ierr); 1366606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 1367606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 1368d5b485abSSatish Balay 1369606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 1370606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 1371d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1372606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1373d5b485abSSatish Balay } 1374d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1375606d414cSSatish Balay ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1376606d414cSSatish Balay ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 1377d5b485abSSatish Balay } 1378d5b485abSSatish Balay 1379606d414cSSatish Balay ierr = PetscFree(sbuf2);CHKERRQ(ierr); 1380606d414cSSatish Balay ierr = PetscFree(rbuf3);CHKERRQ(ierr); 1381606d414cSSatish Balay ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1382606d414cSSatish Balay ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1383606d414cSSatish Balay ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1384606d414cSSatish Balay ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1385606d414cSSatish Balay ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1386d5b485abSSatish Balay 1387233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1388ff0794f7SSatish Balay for (i=0; i<ismax; i++){ 1389233c2ff6SSatish Balay ierr = PetscTableDelete(rowmaps[i]);CHKERRQ(ierr); 1390233c2ff6SSatish Balay ierr = PetscTableDelete(colmaps[i]);CHKERRQ(ierr); 1391233c2ff6SSatish Balay } 1392233c2ff6SSatish Balay ierr = PetscFree(colmaps);CHKERRQ(ierr); 1393233c2ff6SSatish Balay ierr = PetscFree(rowmaps);CHKERRQ(ierr); 1394233c2ff6SSatish Balay #else 1395606d414cSSatish Balay ierr = PetscFree(rmap);CHKERRQ(ierr); 1396233c2ff6SSatish Balay ierr = PetscFree(cmap);CHKERRQ(ierr); 1397233c2ff6SSatish Balay #endif 1398606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1399d5b485abSSatish Balay 1400d5b485abSSatish Balay for (i=0; i<ismax; i++) { 140136f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 140236f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1403d5b485abSSatish Balay } 140434e6ae18SSatish Balay 14053a40ed3dSBarry Smith PetscFunctionReturn(0); 1406d5b485abSSatish Balay } 1407ca161407SBarry Smith 1408