123bdbc58SBarry Smith 2d5b485abSSatish Balay /* 3d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 4d5b485abSSatish Balay and to find submatrices that were shared across processors. 5d5b485abSSatish Balay */ 6c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h> 7c6db04a5SJed Brown #include <petscbt.h> 8d5b485abSSatish Balay 9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char**,PetscInt*,PetscInt**); 10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt**,PetscInt**,PetscInt*); 1109573ac7SBarry Smith extern PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 1209573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 13d5b485abSSatish Balay 144a2ae208SSatish Balay #undef __FUNCT__ 154a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ" 16b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov) 17d5b485abSSatish Balay { 186849ba73SBarry Smith PetscErrorCode ierr; 19d0f46423SBarry Smith PetscInt i,N=C->cmap->N, bs=C->rmap->bs; 2036f4e84dSSatish Balay IS *is_new; 2136f4e84dSSatish Balay 223a40ed3dSBarry Smith PetscFunctionBegin; 23785e854fSJed Brown ierr = PetscMalloc1(imax,&is_new);CHKERRQ(ierr); 24df36731bSBarry Smith /* Convert the indices into block format */ 2505d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr); 2626fbe8dcSKarl Rupp if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n"); 27d5b485abSSatish Balay for (i=0; i<ov; ++i) { 2836f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 29d5b485abSSatish Balay } 306bf464f9SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 3105d8c843SHong Zhang ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr); 326bf464f9SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 33606d414cSSatish Balay ierr = PetscFree(is_new);CHKERRQ(ierr); 343a40ed3dSBarry Smith PetscFunctionReturn(0); 35d5b485abSSatish Balay } 36d5b485abSSatish Balay 37d5b485abSSatish Balay /* 38d5b485abSSatish Balay Sample message format: 39d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 400e9b0e7eSHong Zhang to index sets is[1], is[5] 41d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 42d5b485abSSatish Balay ----------- 43d5b485abSSatish Balay mesg [1] = 1 => is[1] 44d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 45d5b485abSSatish Balay ----------- 46d5b485abSSatish Balay mesg [5] = 5 => is[5] 47d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 48d5b485abSSatish Balay ----------- 49d5b485abSSatish Balay mesg [7] 500e9b0e7eSHong Zhang mesg [n] data(is[1]) 51d5b485abSSatish Balay ----------- 52d5b485abSSatish Balay mesg[n+1] 53d5b485abSSatish Balay mesg[m] data(is[5]) 54d5b485abSSatish Balay ----------- 55d5b485abSSatish Balay 56d5b485abSSatish Balay Notes: 57d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 58d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 59d5b485abSSatish Balay */ 604a2ae208SSatish Balay #undef __FUNCT__ 614a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once" 62db41cccfSHong Zhang PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[]) 63d5b485abSSatish Balay { 64df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 655d0c19d7SBarry Smith const PetscInt **idx,*idx_i; 6624c049a4SHong Zhang PetscInt *n,*w3,*w4,**data,len; 676849ba73SBarry Smith PetscErrorCode ierr; 68b24ad042SBarry Smith PetscMPIInt size,rank,tag1,tag2,*w2,*w1,nrqr; 69245d216aSHong Zhang PetscInt Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr; 708f9f447aSHong Zhang PetscInt *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p; 71b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1,*onodes2,*olengths2; 72f1af5d2fSBarry Smith PetscBT *table; 73d5b485abSSatish Balay MPI_Comm comm; 74d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 75d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 768f9f447aSHong Zhang char *t_p; 77d5b485abSSatish Balay 783a40ed3dSBarry Smith PetscFunctionBegin; 79ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 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 87dcca6d9dSJed Brown ierr = PetscMalloc2(imax+1,&idx,imax,&n);CHKERRQ(ierr); 8824c049a4SHong Zhang 89d5b485abSSatish Balay for (i=0; i<imax; i++) { 90d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 91b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 92d5b485abSSatish Balay } 93d5b485abSSatish Balay 94d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 95d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 96884e879aSBarry Smith ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 97d5b485abSSatish Balay for (i=0; i<imax; i++) { 98b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 99d5b485abSSatish Balay idx_i = idx[i]; 100d5b485abSSatish Balay len = n[i]; 101d5b485abSSatish Balay for (j=0; j<len; j++) { 102d5b485abSSatish Balay row = idx_i[j]; 103f23aa3ddSBarry Smith if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 104ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 105d5b485abSSatish Balay w4[proc]++; 106d5b485abSSatish Balay } 107d5b485abSSatish Balay for (j=0; j<size; j++) { 108d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 109d5b485abSSatish Balay } 110d5b485abSSatish Balay } 111d5b485abSSatish Balay 112d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 113d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1140e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 115d5b485abSSatish Balay w3[rank] = 0; 116d5b485abSSatish Balay for (i=0; i<size; i++) { 117d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 118d5b485abSSatish Balay } 119d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 120854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); 121d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 122d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 123d5b485abSSatish Balay } 124d5b485abSSatish Balay 125d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 126d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 127d5b485abSSatish Balay j = pa[i]; 128d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 129d5b485abSSatish Balay msz += w1[j]; 130d5b485abSSatish Balay } 131d5b485abSSatish Balay 132c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 133a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 134c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 135d5b485abSSatish Balay 136c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 137c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 138d5b485abSSatish Balay 139d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 140dcca6d9dSJed Brown ierr = PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr);CHKERRQ(ierr); 14123bdbc58SBarry Smith ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr); 14223bdbc58SBarry Smith ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 143d5b485abSSatish Balay { 144b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 145d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 146d5b485abSSatish Balay j = pa[i]; 147d5b485abSSatish Balay iptr += ict; 148d5b485abSSatish Balay outdat[j] = iptr; 149d5b485abSSatish Balay ict = w1[j]; 150d5b485abSSatish Balay } 151d5b485abSSatish Balay } 152d5b485abSSatish Balay 153d5b485abSSatish Balay /* Form the outgoing messages */ 154d5b485abSSatish Balay /*plug in the headers*/ 155d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 156d5b485abSSatish Balay j = pa[i]; 157d5b485abSSatish Balay outdat[j][0] = 0; 158b24ad042SBarry Smith ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 159d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 160d5b485abSSatish Balay } 161d5b485abSSatish Balay 162d5b485abSSatish Balay /* Memory for doing local proc's work*/ 163d5b485abSSatish Balay { 1641795a4d1SJed Brown ierr = PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mbs*imax,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,&t_p);CHKERRQ(ierr); 165d5b485abSSatish Balay 166bbd702dbSSatish Balay for (i=0; i<imax; i++) { 167b6410449SSatish Balay table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 168bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 169d5b485abSSatish Balay } 170d5b485abSSatish Balay } 171d5b485abSSatish Balay 172d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 173d5b485abSSatish Balay { 174b24ad042SBarry Smith PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 175f1af5d2fSBarry Smith PetscBT table_i; 176d5b485abSSatish Balay 177d5b485abSSatish Balay for (i=0; i<imax; i++) { 178b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 179d5b485abSSatish Balay n_i = n[i]; 180d5b485abSSatish Balay table_i = table[i]; 181d5b485abSSatish Balay idx_i = idx[i]; 182d5b485abSSatish Balay data_i = data[i]; 183d5b485abSSatish Balay isz_i = isz[i]; 184d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 185d5b485abSSatish Balay row = idx_i[j]; 186ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 187d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 188d5b485abSSatish Balay ctr[proc]++; 189d5b485abSSatish Balay *ptr[proc] = row; 190d5b485abSSatish Balay ptr[proc]++; 191d6b45a43SBarry Smith } else { /* Update the local table */ 19226fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 193d5b485abSSatish Balay } 194d5b485abSSatish Balay } 195d5b485abSSatish Balay /* Update the headers for the current IS */ 196d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 197d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 198d5b485abSSatish Balay outdat_j = outdat[j]; 199d5b485abSSatish Balay k = ++outdat_j[0]; 200d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 201d5b485abSSatish Balay outdat_j[2*k-1] = i; 202d5b485abSSatish Balay } 203d5b485abSSatish Balay } 204d5b485abSSatish Balay isz[i] = isz_i; 205d5b485abSSatish Balay } 206d5b485abSSatish Balay } 207d5b485abSSatish Balay 208d5b485abSSatish Balay /* Now post the sends */ 209854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 210d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 211d5b485abSSatish Balay j = pa[i]; 212b24ad042SBarry Smith ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 213d5b485abSSatish Balay } 214d5b485abSSatish Balay 215d5b485abSSatish Balay /* No longer need the original indices*/ 216d5b485abSSatish Balay for (i=0; i<imax; ++i) { 217d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 218d5b485abSSatish Balay } 21924c049a4SHong Zhang ierr = PetscFree2(idx,n);CHKERRQ(ierr); 220d5b485abSSatish Balay 221d5b485abSSatish Balay for (i=0; i<imax; ++i) { 2226bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 223d5b485abSSatish Balay } 224d5b485abSSatish Balay 225d5b485abSSatish Balay /* Do Local work*/ 226df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 227d5b485abSSatish Balay 228d5b485abSSatish Balay /* Receive messages*/ 229854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&recv_status);CHKERRQ(ierr); 2300c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);} 231d5b485abSSatish Balay 232854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_status);CHKERRQ(ierr); 2330c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 234d5b485abSSatish Balay 235d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 23623bdbc58SBarry Smith ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr); 23723bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 238d5b485abSSatish Balay 239854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&xdata);CHKERRQ(ierr); 240854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&isz1);CHKERRQ(ierr); 241df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 242c05d87d6SBarry Smith ierr = PetscFree(rbuf[0]);CHKERRQ(ierr); 243606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 244d5b485abSSatish Balay 245d5b485abSSatish Balay /* Send the data back*/ 246d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 247d5b485abSSatish Balay { 248b24ad042SBarry Smith PetscMPIInt *rw1; 249d5b485abSSatish Balay 250884e879aSBarry Smith ierr = PetscCalloc1(size,&rw1);CHKERRQ(ierr); 251c7dd2924SSatish Balay 252d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 253d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 254e32f2f54SBarry Smith if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 255d5b485abSSatish Balay rw1[proc] = isz1[i]; 256d5b485abSSatish Balay } 257d5b485abSSatish Balay 258c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 259c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 260d5b485abSSatish Balay 261c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 262c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 26303512dc5SSatish Balay ierr = PetscFree(rw1);CHKERRQ(ierr); 264c7dd2924SSatish Balay } 265c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 266c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 267d5b485abSSatish Balay 268d5b485abSSatish Balay /* Now post the sends */ 269854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 270d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 271d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 272b24ad042SBarry Smith ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 273d5b485abSSatish Balay } 274d5b485abSSatish Balay 275d5b485abSSatish Balay /* receive work done on other processors*/ 276d5b485abSSatish Balay { 277b24ad042SBarry Smith PetscMPIInt idex; 278b24ad042SBarry Smith PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 279f1af5d2fSBarry Smith PetscBT table_i; 280d5b485abSSatish Balay MPI_Status *status2; 281d5b485abSSatish Balay 282854ce69bSBarry Smith ierr = PetscMalloc1(PetscMax(nrqr,nrqs)+1,&status2);CHKERRQ(ierr); 283d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 284999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 285d5b485abSSatish Balay /* Process the message*/ 286999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 287d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 288999d9058SBarry Smith jmax = rbuf2[idex][0]; 289d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 290d5b485abSSatish Balay max = rbuf2_i[2*j]; 291d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 292d5b485abSSatish Balay isz_i = isz[is_no]; 293d5b485abSSatish Balay data_i = data[is_no]; 294d5b485abSSatish Balay table_i = table[is_no]; 295d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 296d5b485abSSatish Balay row = rbuf2_i[ct1]; 29726fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 298d5b485abSSatish Balay } 299d5b485abSSatish Balay isz[is_no] = isz_i; 300d5b485abSSatish Balay } 301d5b485abSSatish Balay } 3020c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);} 303606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 304d5b485abSSatish Balay } 305d5b485abSSatish Balay 306d5b485abSSatish Balay for (i=0; i<imax; ++i) { 30770b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 308d5b485abSSatish Balay } 309d5b485abSSatish Balay 310c7dd2924SSatish Balay 311c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 312c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 313c7dd2924SSatish Balay 314606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 315c05d87d6SBarry Smith ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr); 316606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 317606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 318606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 319606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 320606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 3218f9f447aSHong Zhang ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr); 322606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 323606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 324606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 325606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 326606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 3273a40ed3dSBarry Smith PetscFunctionReturn(0); 328d5b485abSSatish Balay } 329d5b485abSSatish Balay 3304a2ae208SSatish Balay #undef __FUNCT__ 3314a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local" 332d5b485abSSatish Balay /* 333df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 334d5b485abSSatish Balay the work on the local processor. 335d5b485abSSatish Balay 336d5b485abSSatish Balay Inputs: 337df36731bSBarry Smith C - MAT_MPIBAIJ; 338d5b485abSSatish Balay imax - total no of index sets processed at a time; 339df36731bSBarry Smith table - an array of char - size = Mbs bits. 340d5b485abSSatish Balay 341d5b485abSSatish Balay Output: 3420e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 343d5b485abSSatish Balay to each index set; 344d5b485abSSatish Balay data - pointer to the solutions 345d5b485abSSatish Balay */ 346b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 347d5b485abSSatish Balay { 348df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 349d5b485abSSatish Balay Mat A = c->A,B = c->B; 350df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 351b24ad042SBarry Smith PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 352b24ad042SBarry Smith PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 353f1af5d2fSBarry Smith PetscBT table_i; 354d5b485abSSatish Balay 3553a40ed3dSBarry Smith PetscFunctionBegin; 356899cda47SBarry Smith rstart = c->rstartbs; 357899cda47SBarry Smith cstart = c->cstartbs; 358d5b485abSSatish Balay ai = a->i; 359df36731bSBarry Smith aj = a->j; 360d5b485abSSatish Balay bi = b->i; 361df36731bSBarry Smith bj = b->j; 362d5b485abSSatish Balay garray = c->garray; 363d5b485abSSatish Balay 364d5b485abSSatish Balay 365d5b485abSSatish Balay for (i=0; i<imax; i++) { 366d5b485abSSatish Balay data_i = data[i]; 367d5b485abSSatish Balay table_i = table[i]; 368d5b485abSSatish Balay isz_i = isz[i]; 369d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 370d5b485abSSatish Balay row = data_i[j] - rstart; 371d5b485abSSatish Balay start = ai[row]; 372d5b485abSSatish Balay end = ai[row+1]; 373d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 374df36731bSBarry Smith val = aj[k] + cstart; 37526fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 376d5b485abSSatish Balay } 377d5b485abSSatish Balay start = bi[row]; 378d5b485abSSatish Balay end = bi[row+1]; 379d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 380df36731bSBarry Smith val = garray[bj[k]]; 38126fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 382d5b485abSSatish Balay } 383d5b485abSSatish Balay } 384d5b485abSSatish Balay isz[i] = isz_i; 385d5b485abSSatish Balay } 3863a40ed3dSBarry Smith PetscFunctionReturn(0); 387d5b485abSSatish Balay } 3884a2ae208SSatish Balay #undef __FUNCT__ 3894a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive" 390d5b485abSSatish Balay /* 391df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 392d5b485abSSatish Balay and return the output 393d5b485abSSatish Balay 394d5b485abSSatish Balay Input: 395d5b485abSSatish Balay C - the matrix 396d5b485abSSatish Balay nrqr - no of messages being processed. 397d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 398d5b485abSSatish Balay 399d5b485abSSatish Balay Output: 400d5b485abSSatish Balay xdata - array of messages to be sent back 401d5b485abSSatish Balay isz1 - size of each message 402d5b485abSSatish Balay 403a8c7a070SBarry Smith For better efficiency perhaps we should malloc separately each xdata[i], 404d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 4050e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of 406d5b485abSSatish Balay memory is used. 407d5b485abSSatish Balay 408d5b485abSSatish Balay */ 409b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 410d5b485abSSatish Balay { 411df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 412d5b485abSSatish Balay Mat A = c->A,B = c->B; 413df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 4146849ba73SBarry Smith PetscErrorCode ierr; 415b24ad042SBarry Smith PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 416b24ad042SBarry Smith PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 417d2a212eaSBarry Smith PetscInt val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 418b24ad042SBarry Smith PetscInt *rbuf_i,kmax,rbuf_0; 419f1af5d2fSBarry Smith PetscBT xtable; 420d5b485abSSatish Balay 4213a40ed3dSBarry Smith PetscFunctionBegin; 422df36731bSBarry Smith Mbs = c->Mbs; 423899cda47SBarry Smith rstart = c->rstartbs; 424899cda47SBarry Smith cstart = c->cstartbs; 425d5b485abSSatish Balay ai = a->i; 426df36731bSBarry Smith aj = a->j; 427d5b485abSSatish Balay bi = b->i; 428df36731bSBarry Smith bj = b->j; 429d5b485abSSatish Balay garray = c->garray; 430d5b485abSSatish Balay 431d5b485abSSatish Balay 432d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 433d5b485abSSatish Balay rbuf_i = rbuf[i]; 434d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 435d5b485abSSatish Balay ct += rbuf_0; 43626fbe8dcSKarl Rupp for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j]; 437d5b485abSSatish Balay } 438d5b485abSSatish Balay 439701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 440701b8814SBarry Smith else max1 = 1; 441d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 442785e854fSJed Brown ierr = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr); 443d5b485abSSatish Balay ++no_malloc; 44453b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr); 445b24ad042SBarry Smith ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr); 446d5b485abSSatish Balay 447d5b485abSSatish Balay ct3 = 0; 448d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 449d5b485abSSatish Balay rbuf_i = rbuf[i]; 450d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 451d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 452d5b485abSSatish Balay ct2 = ct1; 453d5b485abSSatish Balay ct3 += ct1; 454d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 4556831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 456d5b485abSSatish Balay oct2 = ct2; 457d5b485abSSatish Balay kmax = rbuf_i[2*j]; 458d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 459d5b485abSSatish Balay row = rbuf_i[ct1]; 460f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 461d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 462b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 463785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 464b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 465606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 466d5b485abSSatish Balay xdata[0] = tmp; 467d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 46826fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 469d5b485abSSatish Balay } 470d5b485abSSatish Balay xdata[i][ct2++] = row; 471d5b485abSSatish Balay ct3++; 472d5b485abSSatish Balay } 473d5b485abSSatish Balay } 474d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 475d5b485abSSatish Balay row = xdata[i][k] - rstart; 476d5b485abSSatish Balay start = ai[row]; 477d5b485abSSatish Balay end = ai[row+1]; 478d5b485abSSatish Balay for (l=start; l<end; l++) { 479df36731bSBarry Smith val = aj[l] + cstart; 480f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 481d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 482b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 483785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 484b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 485606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 486d5b485abSSatish Balay xdata[0] = tmp; 487d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 48826fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 489d5b485abSSatish Balay } 490d5b485abSSatish Balay xdata[i][ct2++] = val; 491d5b485abSSatish Balay ct3++; 492d5b485abSSatish Balay } 493d5b485abSSatish Balay } 494d5b485abSSatish Balay start = bi[row]; 495d5b485abSSatish Balay end = bi[row+1]; 496d5b485abSSatish Balay for (l=start; l<end; l++) { 497df36731bSBarry Smith val = garray[bj[l]]; 498f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 499d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 500b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 501785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 502b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 503606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 504d5b485abSSatish Balay xdata[0] = tmp; 505d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 50626fbe8dcSKarl Rupp for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 507d5b485abSSatish Balay } 508d5b485abSSatish Balay xdata[i][ct2++] = val; 509d5b485abSSatish Balay ct3++; 510d5b485abSSatish Balay } 511d5b485abSSatish Balay } 512d5b485abSSatish Balay } 513d5b485abSSatish Balay /* Update the header*/ 514d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 515d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 516d5b485abSSatish Balay } 517d5b485abSSatish Balay xdata[i][0] = rbuf_0; 518d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 519d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 520d5b485abSSatish Balay } 52194bacf5dSBarry Smith ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr); 5221e2582c4SBarry Smith ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr); 5233a40ed3dSBarry Smith PetscFunctionReturn(0); 524d5b485abSSatish Balay } 525d5b485abSSatish Balay 5264a2ae208SSatish Balay #undef __FUNCT__ 5274a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ" 528b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 529d5b485abSSatish Balay { 53036f4e84dSSatish Balay IS *isrow_new,*iscol_new; 531cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 5326849ba73SBarry Smith PetscErrorCode ierr; 5334da72fa9SHong Zhang PetscInt nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs; 5344da72fa9SHong Zhang PetscBool colflag,*allcolumns,*allrows; 535a2feddc5SSatish Balay 5363a40ed3dSBarry Smith PetscFunctionBegin; 53729dcf524SDmitry Karpeev /* The compression and expansion should be avoided. Doesn't point 53829dcf524SDmitry Karpeev out errors, might change the indices, hence buggey */ 539dcca6d9dSJed Brown ierr = PetscMalloc2(ismax+1,&isrow_new,ismax+1,&iscol_new);CHKERRQ(ierr); 54005d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_new);CHKERRQ(ierr); 54105d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_new);CHKERRQ(ierr); 542cf4f527aSSatish Balay 543307b7a18SHong Zhang /* Check for special case: each processor gets entire matrix columns */ 544dcca6d9dSJed Brown ierr = PetscMalloc2(ismax+1,&allcolumns,ismax+1,&allrows);CHKERRQ(ierr); 545307b7a18SHong Zhang for (i=0; i<ismax; i++) { 546307b7a18SHong Zhang ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 547307b7a18SHong Zhang ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr); 54826fbe8dcSKarl Rupp if (colflag && ncol == C->cmap->N) allcolumns[i] = PETSC_TRUE; 54926fbe8dcSKarl Rupp else allcolumns[i] = PETSC_FALSE; 5504da72fa9SHong Zhang 5514da72fa9SHong Zhang ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 5524da72fa9SHong Zhang ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr); 55326fbe8dcSKarl Rupp if (colflag && nrow == C->rmap->N) allrows[i] = PETSC_TRUE; 55426fbe8dcSKarl Rupp else allrows[i] = PETSC_FALSE; 555307b7a18SHong Zhang } 556307b7a18SHong Zhang 557cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 558cf4f527aSSatish Balay if (scall != MAT_REUSE_MATRIX) { 559854ce69bSBarry Smith ierr = PetscMalloc1(ismax+1,submat);CHKERRQ(ierr); 560cf4f527aSSatish Balay } 561cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 562b24ad042SBarry Smith nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 563329f5518SBarry Smith if (!nmax) nmax = 1; 564cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); 565cf4f527aSSatish Balay 566653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 567*b2566f29SBarry Smith ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 568cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 569cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 570cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 571cf4f527aSSatish Balay else max_no = ismax-pos; 5724da72fa9SHong Zhang ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 573cf4f527aSSatish Balay pos += max_no; 574cf4f527aSSatish Balay } 57536f4e84dSSatish Balay 57636f4e84dSSatish Balay for (i=0; i<ismax; i++) { 5776bf464f9SBarry Smith ierr = ISDestroy(&isrow_new[i]);CHKERRQ(ierr); 5786bf464f9SBarry Smith ierr = ISDestroy(&iscol_new[i]);CHKERRQ(ierr); 57936f4e84dSSatish Balay } 58023bdbc58SBarry Smith ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr); 5814da72fa9SHong Zhang ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr); 5823a40ed3dSBarry Smith PetscFunctionReturn(0); 583a2feddc5SSatish Balay } 584a2feddc5SSatish Balay 585233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 5864a2ae208SSatish Balay #undef __FUNCT__ 5874a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc" 588e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 589233c2ff6SSatish Balay { 590e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 5914dc2109aSBarry Smith PetscMPIInt fproc; 5924dc2109aSBarry Smith PetscErrorCode ierr; 593233c2ff6SSatish Balay 594233c2ff6SSatish Balay PetscFunctionBegin; 5954dc2109aSBarry Smith ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr); 59623ce1328SBarry Smith if (fproc > size) fproc = size; 597e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 598e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 599233c2ff6SSatish Balay else fproc++; 600233c2ff6SSatish Balay } 601e005ede5SBarry Smith *rank = fproc; 602233c2ff6SSatish Balay PetscFunctionReturn(0); 603233c2ff6SSatish Balay } 604233c2ff6SSatish Balay #endif 605233c2ff6SSatish Balay 606a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 607b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 6084a2ae208SSatish Balay #undef __FUNCT__ 6094a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local" 6104da72fa9SHong Zhang PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 611a2feddc5SSatish Balay { 612df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 613cf4f527aSSatish Balay Mat A = c->A; 614df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 6155d0c19d7SBarry Smith const PetscInt **irow,**icol,*irow_i; 6165d0c19d7SBarry Smith PetscInt *nrow,*ncol,*w3,*w4,start; 6176849ba73SBarry Smith PetscErrorCode ierr; 6189405f653SBarry Smith PetscMPIInt size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc; 6199405f653SBarry Smith PetscInt **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row; 620b24ad042SBarry Smith PetscInt nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol; 621b24ad042SBarry Smith PetscInt **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 6225d0c19d7SBarry Smith PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 623052f0c41SBarry Smith PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 624d0f46423SBarry Smith PetscInt bs =C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 625899cda47SBarry Smith PetscInt cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 626899cda47SBarry Smith PetscInt *bmap = c->garray,ctmp,rstart=c->rstartbs; 6277a868f3eSHong Zhang MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3; 6287a868f3eSHong Zhang MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3; 629d5b485abSSatish Balay MPI_Comm comm; 630ace3abfcSBarry Smith PetscBool flag; 631b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1; 6327a868f3eSHong Zhang PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 63326fbe8dcSKarl Rupp 6347a868f3eSHong Zhang /* variables below are used for the matrix numerical values - case of !ijonly */ 6357a868f3eSHong Zhang MPI_Request *r_waits4,*s_waits4; 6367a868f3eSHong Zhang MPI_Status *r_status4,*s_status4; 6370298fd71SBarry Smith MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a = NULL,*sbuf_aa_i,*vworkA = NULL,*vworkB = NULL; 6387a868f3eSHong Zhang MatScalar *a_a=a->a,*b_a=b->a; 639c7dd2924SSatish Balay 64080d608c0SSatish Balay #if defined(PETSC_USE_CTABLE) 641b24ad042SBarry Smith PetscInt tt; 6420298fd71SBarry Smith PetscTable *rmap,*cmap,rmap_i,cmap_i=NULL; 643233c2ff6SSatish Balay #else 6440298fd71SBarry Smith PetscInt **cmap,*cmap_i=NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 645233c2ff6SSatish Balay #endif 646d5b485abSSatish Balay 6473a40ed3dSBarry Smith PetscFunctionBegin; 648ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 6497adad957SLisandro Dalcin tag0 = ((PetscObject)C)->tag; 650d5b485abSSatish Balay size = c->size; 651d5b485abSSatish Balay rank = c->rank; 652d5b485abSSatish Balay 65334e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 65434e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 65534e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 65634e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 65734e6ae18SSatish Balay 658052f0c41SBarry Smith #if defined(PETSC_USE_CTABLE) 659dcca6d9dSJed Brown ierr = PetscMalloc4(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol);CHKERRQ(ierr); 660052f0c41SBarry Smith #else 661dcca6d9dSJed Brown ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,Mbs+1,&rtable);CHKERRQ(ierr); 662d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 663d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 66482a7e548SBarry Smith jmax = C->rmap->range[i+1]/bs; 66526fbe8dcSKarl Rupp for (; j<jmax; j++) rtable[j] = i; 666d5b485abSSatish Balay } 667233c2ff6SSatish Balay #endif 668233c2ff6SSatish Balay 669233c2ff6SSatish Balay for (i=0; i<ismax; i++) { 6704da72fa9SHong Zhang if (allrows[i]) { 6710298fd71SBarry Smith irow[i] = NULL; 6724da72fa9SHong Zhang nrow[i] = C->rmap->N/bs; 6734da72fa9SHong Zhang } else { 674233c2ff6SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 675233c2ff6SSatish Balay ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 6764da72fa9SHong Zhang } 6774da72fa9SHong Zhang 678307b7a18SHong Zhang if (allcolumns[i]) { 6790298fd71SBarry Smith icol[i] = NULL; 680307b7a18SHong Zhang ncol[i] = C->cmap->N/bs; 681307b7a18SHong Zhang } else { 682307b7a18SHong Zhang ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 683233c2ff6SSatish Balay ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 684233c2ff6SSatish Balay } 685307b7a18SHong Zhang } 686d5b485abSSatish Balay 687d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg,and buffer space 688d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 689884e879aSBarry Smith ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 690d5b485abSSatish Balay for (i=0; i<ismax; i++) { 691b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 692d5b485abSSatish Balay jmax = nrow[i]; 693d5b485abSSatish Balay irow_i = irow[i]; 694d5b485abSSatish Balay for (j=0; j<jmax; j++) { 69526fbe8dcSKarl Rupp if (allrows[i]) row = j; 69626fbe8dcSKarl Rupp else row = irow_i[j]; 69726fbe8dcSKarl Rupp 698233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 699899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 700233c2ff6SSatish Balay #else 701d5b485abSSatish Balay proc = rtable[row]; 702233c2ff6SSatish Balay #endif 703d5b485abSSatish Balay w4[proc]++; 704d5b485abSSatish Balay } 705d5b485abSSatish Balay for (j=0; j<size; j++) { 706d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 707d5b485abSSatish Balay } 708d5b485abSSatish Balay } 709d5b485abSSatish Balay 710d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 711e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 712d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 713d5b485abSSatish Balay w3[rank] = 0; 714d5b485abSSatish Balay for (i=0; i<size; i++) { 715d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 716d5b485abSSatish Balay } 717854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 718d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 719d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 720d5b485abSSatish Balay } 721d5b485abSSatish Balay 722d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 723d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 724d5b485abSSatish Balay j = pa[i]; 725d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 726d5b485abSSatish Balay msz += w1[j]; 727d5b485abSSatish Balay } 728d5b485abSSatish Balay 729c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 730a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 731c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 732d5b485abSSatish Balay 733c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 734c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 735c7dd2924SSatish Balay 736c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 737c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 738d5b485abSSatish Balay 739d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 740dcca6d9dSJed Brown ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 74123bdbc58SBarry Smith ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 74223bdbc58SBarry Smith ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 743d5b485abSSatish Balay { 744b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 745d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 746d5b485abSSatish Balay j = pa[i]; 747d5b485abSSatish Balay iptr += ict; 748d5b485abSSatish Balay sbuf1[j] = iptr; 749d5b485abSSatish Balay ict = w1[j]; 750d5b485abSSatish Balay } 751d5b485abSSatish Balay } 752d5b485abSSatish Balay 753d5b485abSSatish Balay /* Form the outgoing messages */ 754d5b485abSSatish Balay /* Initialise the header space */ 755d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 756d5b485abSSatish Balay j = pa[i]; 757d5b485abSSatish Balay sbuf1[j][0] = 0; 758b24ad042SBarry Smith ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 759d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 760d5b485abSSatish Balay } 761d5b485abSSatish Balay 762d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 763d5b485abSSatish Balay for (i=0; i<ismax; i++) { 764b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 765d5b485abSSatish Balay irow_i = irow[i]; 766d5b485abSSatish Balay jmax = nrow[i]; 767d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 76826fbe8dcSKarl Rupp if (allrows[i]) row = j; 76926fbe8dcSKarl Rupp else row = irow_i[j]; 77026fbe8dcSKarl Rupp 771233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 772899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 773233c2ff6SSatish Balay #else 774d5b485abSSatish Balay proc = rtable[row]; 775233c2ff6SSatish Balay #endif 776d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 777d5b485abSSatish Balay ctr[proc]++; 778d5b485abSSatish Balay *ptr[proc] = row; 779d5b485abSSatish Balay ptr[proc]++; 780d5b485abSSatish Balay } 781d5b485abSSatish Balay } 782d5b485abSSatish Balay /* Update the headers for the current IS */ 783d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 78406ef35abSSatish Balay if ((ctr_j = ctr[j])) { 785d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 786d5b485abSSatish Balay k = ++sbuf1_j[0]; 787d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 788d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 789d5b485abSSatish Balay } 790d5b485abSSatish Balay } 791d5b485abSSatish Balay } 792d5b485abSSatish Balay 793d5b485abSSatish Balay /* Now post the sends */ 794854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 795d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 796d5b485abSSatish Balay j = pa[i]; 797b24ad042SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 798d5b485abSSatish Balay } 799d5b485abSSatish Balay 800d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 801854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 802854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&rbuf2);CHKERRQ(ierr); 803d5b485abSSatish Balay rbuf2[0] = tmp + msz; 804d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 805d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 806d5b485abSSatish Balay } 807d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 808d5b485abSSatish Balay j = pa[i]; 809b24ad042SBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 810d5b485abSSatish Balay } 811d5b485abSSatish Balay 812d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 813d5b485abSSatish Balay 814d5b485abSSatish Balay /* Receive messages*/ 815854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 816854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 817dcca6d9dSJed Brown ierr = PetscMalloc3(nrqr+1,&sbuf2,nrqr,&req_size,nrqr,&req_source);CHKERRQ(ierr); 818d5b485abSSatish Balay { 819df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 820b24ad042SBarry Smith PetscInt *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 821d5b485abSSatish Balay 822d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 823999d9058SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 82426fbe8dcSKarl Rupp 825999d9058SBarry Smith req_size[idex] = 0; 826999d9058SBarry Smith rbuf1_i = rbuf1[idex]; 827d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 828b24ad042SBarry Smith ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 829785e854fSJed Brown ierr = PetscMalloc1(end,&sbuf2[idex]);CHKERRQ(ierr); 830999d9058SBarry Smith sbuf2_i = sbuf2[idex]; 831d5b485abSSatish Balay for (j=start; j<end; j++) { 832d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 833d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 834d5b485abSSatish Balay sbuf2_i[j] = ncols; 835999d9058SBarry Smith req_size[idex] += ncols; 836d5b485abSSatish Balay } 837999d9058SBarry Smith req_source[idex] = r_status1[i].MPI_SOURCE; 838d5b485abSSatish Balay /* form the header */ 839999d9058SBarry Smith sbuf2_i[0] = req_size[idex]; 84026fbe8dcSKarl Rupp for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 841b24ad042SBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 842d5b485abSSatish Balay } 843d5b485abSSatish Balay } 844606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 845606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 846d5b485abSSatish Balay 847d5b485abSSatish Balay /* recv buffer sizes */ 848d5b485abSSatish Balay /* Receive messages*/ 849854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&rbuf3);CHKERRQ(ierr); 850854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 851854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 8527a868f3eSHong Zhang if (!ijonly) { 853854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 854854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 8557a868f3eSHong Zhang } 856d5b485abSSatish Balay 857d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 858999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 859785e854fSJed Brown ierr = PetscMalloc1(rbuf2[idex][0],&rbuf3[idex]);CHKERRQ(ierr); 860b24ad042SBarry Smith ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 8617a868f3eSHong Zhang if (!ijonly) { 862785e854fSJed Brown ierr = PetscMalloc1(rbuf2[idex][0]*bs2,&rbuf4[idex]);CHKERRQ(ierr); 863b24ad042SBarry Smith ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 864d5b485abSSatish Balay } 8657a868f3eSHong Zhang } 866606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 867606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 868d5b485abSSatish Balay 869d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 870854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 871854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 872d5b485abSSatish Balay 8730c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 8740c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 875606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 876606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 877606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 878606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 879d5b485abSSatish Balay 880d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 881854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 882d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 883854ce69bSBarry Smith ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 884d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 885d5b485abSSatish Balay 886854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 887d5b485abSSatish Balay { 888d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 889d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 890d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 891d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 892d5b485abSSatish Balay ct2 = 0; 893d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 894d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 895d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 896e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 89705b2c859SBarry Smith nzA = a_i[row+1] - a_i[row]; 89805b2c859SBarry Smith nzB = b_i[row+1] - b_i[row]; 899d5b485abSSatish Balay ncols = nzA + nzB; 90005b2c859SBarry Smith cworkA = a_j + a_i[row]; 90105b2c859SBarry Smith cworkB = b_j + b_i[row]; 902d5b485abSSatish Balay 903d5b485abSSatish Balay /* load the column indices for this row into cols*/ 904d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 905d5b485abSSatish Balay for (l=0; l<nzB; l++) { 906d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 907d5b485abSSatish Balay else break; 908d5b485abSSatish Balay } 909d5b485abSSatish Balay imark = l; 910d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 911d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 912d5b485abSSatish Balay ct2 += ncols; 913d5b485abSSatish Balay } 914d5b485abSSatish Balay } 915b24ad042SBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 916d5b485abSSatish Balay } 917d5b485abSSatish Balay } 918854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_status3);CHKERRQ(ierr); 919854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_status3);CHKERRQ(ierr); 920d5b485abSSatish Balay 921d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 9227a868f3eSHong Zhang if (!ijonly) { 923854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 924d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 925785e854fSJed Brown ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 926a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 927d5b485abSSatish Balay 928854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 929d5b485abSSatish Balay { 930d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 931d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 932d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 933d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 934d5b485abSSatish Balay ct2 = 0; 935d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 936d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 937d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 938e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 93905b2c859SBarry Smith nzA = a_i[row+1] - a_i[row]; 94005b2c859SBarry Smith nzB = b_i[row+1] - b_i[row]; 941d5b485abSSatish Balay ncols = nzA + nzB; 94205b2c859SBarry Smith cworkB = b_j + b_i[row]; 94305b2c859SBarry Smith vworkA = a_a + a_i[row]*bs2; 94405b2c859SBarry Smith vworkB = b_a + b_i[row]*bs2; 945d5b485abSSatish Balay 946d5b485abSSatish Balay /* load the column values for this row into vals*/ 9475b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 948d5b485abSSatish Balay for (l=0; l<nzB; l++) { 9493a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 9503eda8832SBarry Smith ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 95126fbe8dcSKarl Rupp } else break; 952d5b485abSSatish Balay } 953d5b485abSSatish Balay imark = l; 9543a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 9553eda8832SBarry Smith ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9563a40ed3dSBarry Smith } 9573a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 9583eda8832SBarry Smith ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9593a40ed3dSBarry Smith } 960d5b485abSSatish Balay ct2 += ncols; 961d5b485abSSatish Balay } 962d5b485abSSatish Balay } 9633eda8832SBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 964d5b485abSSatish Balay } 965d5b485abSSatish Balay } 966854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 967854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 9687a868f3eSHong Zhang } 969533163c2SBarry Smith ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 970606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 971d5b485abSSatish Balay 972d5b485abSSatish Balay /* Form the matrix */ 973307b7a18SHong Zhang /* create col map: global col of C -> local col of submatrices */ 974d5b485abSSatish Balay { 9755d0c19d7SBarry Smith const PetscInt *icol_i; 976233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 977854ce69bSBarry Smith ierr = PetscMalloc1(1+ismax,&cmap);CHKERRQ(ierr); 978ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 979307b7a18SHong Zhang if (!allcolumns[i]) { 9808fa52d88SHong Zhang ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 981307b7a18SHong Zhang jmax = ncol[i]; 982307b7a18SHong Zhang icol_i = icol[i]; 9838fa52d88SHong Zhang cmap_i = cmap[i]; 984307b7a18SHong Zhang for (j=0; j<jmax; j++) { 9853861aac3SJed Brown ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 986307b7a18SHong Zhang } 987307b7a18SHong Zhang } else { 9880298fd71SBarry Smith cmap[i] = NULL; 989307b7a18SHong Zhang } 990233c2ff6SSatish Balay } 991233c2ff6SSatish Balay #else 992785e854fSJed Brown ierr = PetscMalloc1(ismax,&cmap);CHKERRQ(ierr); 993d5b485abSSatish Balay for (i=0; i<ismax; i++) { 9948fa52d88SHong Zhang if (!allcolumns[i]) { 995884e879aSBarry Smith ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 996d5b485abSSatish Balay jmax = ncol[i]; 997d5b485abSSatish Balay icol_i = icol[i]; 998d5b485abSSatish Balay cmap_i = cmap[i]; 999d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1000d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1001d5b485abSSatish Balay } 10028fa52d88SHong Zhang } else { /* allcolumns[i] */ 10030298fd71SBarry Smith cmap[i] = NULL; 10048fa52d88SHong Zhang } 1005d5b485abSSatish Balay } 1006307b7a18SHong Zhang #endif 1007d5b485abSSatish Balay } 1008d5b485abSSatish Balay 1009d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 101026fbe8dcSKarl Rupp for (i=0,j=0; i<ismax; i++) j += nrow[i]; 1011884e879aSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr); 1012b24ad042SBarry Smith lens[0] = (PetscInt*)(lens + ismax); 1013b24ad042SBarry Smith ierr = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr); 101426fbe8dcSKarl Rupp for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 1015d5b485abSSatish Balay 1016d5b485abSSatish Balay /* Update lens from local data */ 1017d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1018d5b485abSSatish Balay jmax = nrow[i]; 10198fa52d88SHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 1020d5b485abSSatish Balay irow_i = irow[i]; 1021d5b485abSSatish Balay lens_i = lens[i]; 1022d5b485abSSatish Balay for (j=0; j<jmax; j++) { 102326fbe8dcSKarl Rupp if (allrows[i]) row = j; 102426fbe8dcSKarl Rupp else row = irow_i[j]; 102526fbe8dcSKarl Rupp 1026233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1027899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1028233c2ff6SSatish Balay #else 1029d5b485abSSatish Balay proc = rtable[row]; 1030233c2ff6SSatish Balay #endif 1031d5b485abSSatish Balay if (proc == rank) { 103236f4e84dSSatish Balay /* Get indices from matA and then from matB */ 103336f4e84dSSatish Balay row = row - rstart; 103405b2c859SBarry Smith nzA = a_i[row+1] - a_i[row]; 103505b2c859SBarry Smith nzB = b_i[row+1] - b_i[row]; 103605b2c859SBarry Smith cworkA = a_j + a_i[row]; 103705b2c859SBarry Smith cworkB = b_j + b_i[row]; 1038307b7a18SHong Zhang if (!allcolumns[i]) { 1039233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1040233c2ff6SSatish Balay for (k=0; k<nzA; k++) { 10418fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 104226fbe8dcSKarl Rupp if (tt) lens_i[j]++; 1043233c2ff6SSatish Balay } 1044233c2ff6SSatish Balay for (k=0; k<nzB; k++) { 10458fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 104626fbe8dcSKarl Rupp if (tt) lens_i[j]++; 1047233c2ff6SSatish Balay } 1048307b7a18SHong Zhang 1049233c2ff6SSatish Balay #else 1050ca161407SBarry Smith for (k=0; k<nzA; k++) { 105126fbe8dcSKarl Rupp if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1052ca161407SBarry Smith } 1053ca161407SBarry Smith for (k=0; k<nzB; k++) { 105426fbe8dcSKarl Rupp if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1055d5b485abSSatish Balay } 1056233c2ff6SSatish Balay #endif 1057307b7a18SHong Zhang } else { /* allcolumns */ 1058307b7a18SHong Zhang lens_i[j] = nzA + nzB; 1059307b7a18SHong Zhang } 1060d5b485abSSatish Balay } 1061a2feddc5SSatish Balay } 1062ca161407SBarry Smith } 1063233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1064d5b485abSSatish Balay /* Create row map*/ 1065854ce69bSBarry Smith ierr = PetscMalloc1(1+ismax,&rmap);CHKERRQ(ierr); 1066ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 10678fa52d88SHong Zhang ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1068233c2ff6SSatish Balay } 1069233c2ff6SSatish Balay #else 1070233c2ff6SSatish Balay /* Create row map*/ 1071884e879aSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr); 1072b24ad042SBarry Smith rmap[0] = (PetscInt*)(rmap + ismax); 1073b24ad042SBarry Smith ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr); 107426fbe8dcSKarl Rupp for (i=1; i<ismax; i++) rmap[i] = rmap[i-1] + Mbs; 1075233c2ff6SSatish Balay #endif 1076d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1077d5b485abSSatish Balay irow_i = irow[i]; 1078d5b485abSSatish Balay jmax = nrow[i]; 1079233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 10808fa52d88SHong Zhang rmap_i = rmap[i]; 1081233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 10824da72fa9SHong Zhang if (allrows[i]) { 10833861aac3SJed Brown ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 10844da72fa9SHong Zhang } else { 10853861aac3SJed Brown ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1086233c2ff6SSatish Balay } 10874da72fa9SHong Zhang } 1088233c2ff6SSatish Balay #else 1089233c2ff6SSatish Balay rmap_i = rmap[i]; 1090d5b485abSSatish Balay for (j=0; j<jmax; j++) { 109126fbe8dcSKarl Rupp if (allrows[i]) rmap_i[j] = j; 109226fbe8dcSKarl Rupp else rmap_i[irow_i[j]] = j; 10934da72fa9SHong Zhang } 1094233c2ff6SSatish Balay #endif 1095d5b485abSSatish Balay } 1096d5b485abSSatish Balay 1097d5b485abSSatish Balay /* Update lens from offproc data */ 1098d5b485abSSatish Balay { 1099b24ad042SBarry Smith PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1100b24ad042SBarry Smith PetscMPIInt ii; 1101d5b485abSSatish Balay 1102d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1103b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr); 1104b24ad042SBarry Smith idex = pa[ii]; 1105999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1106d5b485abSSatish Balay jmax = sbuf1_i[0]; 1107d5b485abSSatish Balay ct1 = 2*jmax+1; 1108d5b485abSSatish Balay ct2 = 0; 1109b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1110b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 1111d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1112d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1113d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1114d5b485abSSatish Balay lens_i = lens[is_no]; 11158fa52d88SHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1116d5b485abSSatish Balay rmap_i = rmap[is_no]; 1117d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1118233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 11198fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1120233c2ff6SSatish Balay row--; 112126fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1122233c2ff6SSatish Balay #else 1123d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1124233c2ff6SSatish Balay #endif 1125d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1126d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1127307b7a18SHong Zhang if (!allcolumns[is_no]) { 1128233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 11298fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 113026fbe8dcSKarl Rupp if (tt) lens_i[row]++; 1131233c2ff6SSatish Balay #else 113226fbe8dcSKarl Rupp if (cmap_i[rbuf3_i[ct2]]) lens_i[row]++; 1133233c2ff6SSatish Balay #endif 1134307b7a18SHong Zhang } else { /* allcolumns */ 1135307b7a18SHong Zhang lens_i[row]++; 1136307b7a18SHong Zhang } 1137d5b485abSSatish Balay } 1138d5b485abSSatish Balay } 1139d5b485abSSatish Balay } 1140d5b485abSSatish Balay } 1141d5b485abSSatish Balay } 1142606d414cSSatish Balay ierr = PetscFree(r_status3);CHKERRQ(ierr); 1143606d414cSSatish Balay ierr = PetscFree(r_waits3);CHKERRQ(ierr); 11440c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 1145606d414cSSatish Balay ierr = PetscFree(s_status3);CHKERRQ(ierr); 1146606d414cSSatish Balay ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1147d5b485abSSatish Balay 1148d5b485abSSatish Balay /* Create the submatrices */ 1149d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 11507a868f3eSHong Zhang if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet"); 1151d5b485abSSatish Balay /* 1152d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1153d5b485abSSatish Balay */ 1154d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1155df36731bSBarry Smith mat = (Mat_SeqBAIJ*)(submats[i]->data); 1156e7e72b3dSBarry Smith if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 1157b24ad042SBarry Smith ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 1158e7e72b3dSBarry Smith if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 1159d5b485abSSatish Balay /* Initial matrix as if empty */ 1160b24ad042SBarry Smith ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr); 116126fbe8dcSKarl Rupp 1162d5f3da31SBarry Smith submats[i]->factortype = C->factortype; 1163d5b485abSSatish Balay } 1164ca161407SBarry Smith } else { 11657a868f3eSHong Zhang PetscInt bs_tmp; 116626fbe8dcSKarl Rupp if (ijonly) bs_tmp = 1; 116726fbe8dcSKarl Rupp else bs_tmp = bs; 1168d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1169f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 11707a868f3eSHong Zhang ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr); 11717adad957SLisandro Dalcin ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 11727a868f3eSHong Zhang ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 11737a868f3eSHong Zhang ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 1174d5b485abSSatish Balay } 1175d5b485abSSatish Balay } 1176d5b485abSSatish Balay 1177d5b485abSSatish Balay /* Assemble the matrices */ 1178d5b485abSSatish Balay /* First assemble the local rows */ 1179d5b485abSSatish Balay { 1180b24ad042SBarry Smith PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 11810298fd71SBarry Smith MatScalar *imat_a = NULL; 1182d5b485abSSatish Balay 1183d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1184df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[i]->data; 1185d5b485abSSatish Balay imat_ilen = mat->ilen; 1186d5b485abSSatish Balay imat_j = mat->j; 1187d5b485abSSatish Balay imat_i = mat->i; 11887a868f3eSHong Zhang if (!ijonly) imat_a = mat->a; 11898fa52d88SHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 1190d5b485abSSatish Balay rmap_i = rmap[i]; 1191d5b485abSSatish Balay irow_i = irow[i]; 1192d5b485abSSatish Balay jmax = nrow[i]; 1193d5b485abSSatish Balay for (j=0; j<jmax; j++) { 119426fbe8dcSKarl Rupp if (allrows[i]) row = j; 119526fbe8dcSKarl Rupp else row = irow_i[j]; 119626fbe8dcSKarl Rupp 1197233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1198899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1199233c2ff6SSatish Balay #else 1200d5b485abSSatish Balay proc = rtable[row]; 1201233c2ff6SSatish Balay #endif 1202d5b485abSSatish Balay if (proc == rank) { 120336f4e84dSSatish Balay row = row - rstart; 120436f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 120536f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 120636f4e84dSSatish Balay cworkA = a_j + a_i[row]; 120736f4e84dSSatish Balay cworkB = b_j + b_i[row]; 12087a868f3eSHong Zhang if (!ijonly) { 120936f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 121036f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 12117a868f3eSHong Zhang } 1212233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 12138fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1214233c2ff6SSatish Balay row--; 121526fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1216233c2ff6SSatish Balay #else 121736f4e84dSSatish Balay row = rmap_i[row + rstart]; 1218233c2ff6SSatish Balay #endif 1219df36731bSBarry Smith mat_i = imat_i[row]; 12207a868f3eSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1221d5b485abSSatish Balay mat_j = imat_j + mat_i; 122236f4e84dSSatish Balay ilen_row = imat_ilen[row]; 122336f4e84dSSatish Balay 122436f4e84dSSatish Balay /* load the column indices for this row into cols*/ 1225307b7a18SHong Zhang if (!allcolumns[i]) { 122636f4e84dSSatish Balay for (l=0; l<nzB; l++) { 122736f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 1228233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 12298fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1230233c2ff6SSatish Balay if (tcol) { 1231233c2ff6SSatish Balay #else 123236f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1233233c2ff6SSatish Balay #endif 1234df36731bSBarry Smith *mat_j++ = tcol - 1; 12353eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1236549d3d68SSatish Balay mat_a += bs2; 1237d5b485abSSatish Balay ilen_row++; 1238d5b485abSSatish Balay } 1239ca161407SBarry Smith } else break; 124036f4e84dSSatish Balay } 124136f4e84dSSatish Balay imark = l; 124236f4e84dSSatish Balay for (l=0; l<nzA; l++) { 1243233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 12448fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1245233c2ff6SSatish Balay if (tcol) { 1246233c2ff6SSatish Balay #else 124736f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 1248233c2ff6SSatish Balay #endif 124936f4e84dSSatish Balay *mat_j++ = tcol - 1; 12507a868f3eSHong Zhang if (!ijonly) { 12513eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1252549d3d68SSatish Balay mat_a += bs2; 12537a868f3eSHong Zhang } 125436f4e84dSSatish Balay ilen_row++; 125536f4e84dSSatish Balay } 125636f4e84dSSatish Balay } 125736f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 1258233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 12598fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1260233c2ff6SSatish Balay if (tcol) { 1261233c2ff6SSatish Balay #else 126236f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1263233c2ff6SSatish Balay #endif 126436f4e84dSSatish Balay *mat_j++ = tcol - 1; 12657a868f3eSHong Zhang if (!ijonly) { 12663eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1267549d3d68SSatish Balay mat_a += bs2; 12687a868f3eSHong Zhang } 126936f4e84dSSatish Balay ilen_row++; 127036f4e84dSSatish Balay } 127136f4e84dSSatish Balay } 1272307b7a18SHong Zhang } else { /* allcolumns */ 1273307b7a18SHong Zhang for (l=0; l<nzB; l++) { 1274307b7a18SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1275307b7a18SHong Zhang *mat_j++ = ctmp; 1276307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1277307b7a18SHong Zhang mat_a += bs2; 1278307b7a18SHong Zhang ilen_row++; 1279307b7a18SHong Zhang } else break; 1280307b7a18SHong Zhang } 1281307b7a18SHong Zhang imark = l; 1282307b7a18SHong Zhang for (l=0; l<nzA; l++) { 1283307b7a18SHong Zhang *mat_j++ = cstart+cworkA[l]; 1284307b7a18SHong Zhang if (!ijonly) { 1285307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1286307b7a18SHong Zhang mat_a += bs2; 1287307b7a18SHong Zhang } 1288307b7a18SHong Zhang ilen_row++; 1289307b7a18SHong Zhang } 1290307b7a18SHong Zhang for (l=imark; l<nzB; l++) { 1291307b7a18SHong Zhang *mat_j++ = bmap[cworkB[l]]; 1292307b7a18SHong Zhang if (!ijonly) { 1293307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1294307b7a18SHong Zhang mat_a += bs2; 1295307b7a18SHong Zhang } 1296307b7a18SHong Zhang ilen_row++; 1297307b7a18SHong Zhang } 1298307b7a18SHong Zhang } 1299d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1300d5b485abSSatish Balay } 1301d5b485abSSatish Balay } 1302d5b485abSSatish Balay } 1303d5b485abSSatish Balay } 1304d5b485abSSatish Balay 1305d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1306d5b485abSSatish Balay { 1307b24ad042SBarry Smith PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1308b24ad042SBarry Smith PetscInt *imat_j,*imat_i; 13090298fd71SBarry Smith MatScalar *imat_a = NULL,*rbuf4_i = NULL; 1310b24ad042SBarry Smith PetscMPIInt ii; 1311d5b485abSSatish Balay 1312d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 131326fbe8dcSKarl Rupp if (ijonly) ii = tmp2; 131426fbe8dcSKarl Rupp else { 1315b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr); 13167a868f3eSHong Zhang } 1317b24ad042SBarry Smith idex = pa[ii]; 1318999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1319d5b485abSSatish Balay jmax = sbuf1_i[0]; 1320d5b485abSSatish Balay ct1 = 2*jmax + 1; 1321d5b485abSSatish Balay ct2 = 0; 1322b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1323b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 13247a868f3eSHong Zhang if (!ijonly) rbuf4_i = rbuf4[ii]; 1325d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1326d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 13278fa52d88SHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1328d5b485abSSatish Balay rmap_i = rmap[is_no]; 1329df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[is_no]->data; 1330d5b485abSSatish Balay imat_ilen = mat->ilen; 1331d5b485abSSatish Balay imat_j = mat->j; 1332d5b485abSSatish Balay imat_i = mat->i; 13337a868f3eSHong Zhang if (!ijonly) imat_a = mat->a; 1334d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1335d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1336d5b485abSSatish Balay row = sbuf1_i[ct1]; 1337233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 13388fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 1339233c2ff6SSatish Balay row--; 134026fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1341233c2ff6SSatish Balay #else 1342d5b485abSSatish Balay row = rmap_i[row]; 1343233c2ff6SSatish Balay #endif 1344d5b485abSSatish Balay ilen = imat_ilen[row]; 1345df36731bSBarry Smith mat_i = imat_i[row]; 13467a868f3eSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1347d5b485abSSatish Balay mat_j = imat_j + mat_i; 1348d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1349307b7a18SHong Zhang 1350307b7a18SHong Zhang if (!allcolumns[is_no]) { 1351d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1352233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 13538fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1354233c2ff6SSatish Balay if (tcol) { 1355233c2ff6SSatish Balay #else 1356d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1357233c2ff6SSatish Balay #endif 1358df36731bSBarry Smith *mat_j++ = tcol - 1; 13597a868f3eSHong Zhang if (!ijonly) { 13603eda8832SBarry Smith ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1361549d3d68SSatish Balay mat_a += bs2; 13627a868f3eSHong Zhang } 1363d5b485abSSatish Balay ilen++; 1364d5b485abSSatish Balay } 1365d5b485abSSatish Balay } 1366307b7a18SHong Zhang } else { /* allcolumns */ 1367307b7a18SHong Zhang for (l=0; l<max2; l++,ct2++) { 1368307b7a18SHong Zhang *mat_j++ = rbuf3_i[ct2]; 1369307b7a18SHong Zhang if (!ijonly) { 1370307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1371307b7a18SHong Zhang mat_a += bs2; 1372307b7a18SHong Zhang } 1373307b7a18SHong Zhang ilen++; 1374307b7a18SHong Zhang } 1375307b7a18SHong Zhang } 1376d5b485abSSatish Balay imat_ilen[row] = ilen; 1377d5b485abSSatish Balay } 1378d5b485abSSatish Balay } 1379d5b485abSSatish Balay } 1380d5b485abSSatish Balay } 1381cdc6f3adSToby Isaac /* sort the rows */ 1382cdc6f3adSToby Isaac { 1383cdc6f3adSToby Isaac PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 1384cdc6f3adSToby Isaac MatScalar *imat_a = NULL; 1385cdc6f3adSToby Isaac MatScalar *work; 1386cdc6f3adSToby Isaac 1387cdc6f3adSToby Isaac ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 1388cdc6f3adSToby Isaac for (i=0; i<ismax; i++) { 1389cdc6f3adSToby Isaac mat = (Mat_SeqBAIJ*)submats[i]->data; 1390cdc6f3adSToby Isaac imat_ilen = mat->ilen; 1391cdc6f3adSToby Isaac imat_j = mat->j; 1392cdc6f3adSToby Isaac imat_i = mat->i; 1393cdc6f3adSToby Isaac if (!ijonly) imat_a = mat->a; 1394cdc6f3adSToby Isaac if (allcolumns[i]) continue; 1395cdc6f3adSToby Isaac jmax = nrow[i]; 1396cdc6f3adSToby Isaac for (j=0; j<jmax; j++) { 1397cdc6f3adSToby Isaac mat_i = imat_i[j]; 1398cdc6f3adSToby Isaac if (!ijonly) mat_a = imat_a + mat_i*bs2; 1399cdc6f3adSToby Isaac mat_j = imat_j + mat_i; 1400cdc6f3adSToby Isaac ilen_row = imat_ilen[j]; 1401cdc6f3adSToby Isaac if (!ijonly) {ierr = PetscSortIntWithDataArray(ilen_row,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);} 1402cdc6f3adSToby Isaac else {ierr = PetscSortInt(ilen_row,mat_j);CHKERRQ(ierr);} 1403cdc6f3adSToby Isaac } 1404cdc6f3adSToby Isaac } 1405cdc6f3adSToby Isaac ierr = PetscFree(work);CHKERRQ(ierr); 1406cdc6f3adSToby Isaac } 14077a868f3eSHong Zhang if (!ijonly) { 1408606d414cSSatish Balay ierr = PetscFree(r_status4);CHKERRQ(ierr); 1409606d414cSSatish Balay ierr = PetscFree(r_waits4);CHKERRQ(ierr); 14100c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 1411606d414cSSatish Balay ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1412606d414cSSatish Balay ierr = PetscFree(s_status4);CHKERRQ(ierr); 14137a868f3eSHong Zhang } 1414d5b485abSSatish Balay 1415d5b485abSSatish Balay /* Restore the indices */ 1416d5b485abSSatish Balay for (i=0; i<ismax; i++) { 14174da72fa9SHong Zhang if (!allrows[i]) { 1418d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 14194da72fa9SHong Zhang } 1420307b7a18SHong Zhang if (!allcolumns[i]) { 1421d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1422d5b485abSSatish Balay } 1423307b7a18SHong Zhang } 1424d5b485abSSatish Balay 1425d5b485abSSatish Balay /* Destroy allocated memory */ 142623bdbc58SBarry Smith #if defined(PETSC_USE_CTABLE) 142723bdbc58SBarry Smith ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr); 142823bdbc58SBarry Smith #else 142923bdbc58SBarry Smith ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr); 143023bdbc58SBarry Smith #endif 143123bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 1432606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 1433d5b485abSSatish Balay 143423bdbc58SBarry Smith ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 1435606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 1436606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 1437d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1438606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1439d5b485abSSatish Balay } 1440d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1441606d414cSSatish Balay ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1442d5b485abSSatish Balay } 144323bdbc58SBarry Smith ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr); 1444606d414cSSatish Balay ierr = PetscFree(rbuf3);CHKERRQ(ierr); 1445606d414cSSatish Balay ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1446606d414cSSatish Balay ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 14477a868f3eSHong Zhang if (!ijonly) { 14487a868f3eSHong Zhang for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);} 14497a868f3eSHong Zhang ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1450606d414cSSatish Balay ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1451606d414cSSatish Balay ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 14527a868f3eSHong Zhang } 1453d5b485abSSatish Balay 1454233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1455ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 14568fa52d88SHong Zhang ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr); 1457233c2ff6SSatish Balay } 1458233c2ff6SSatish Balay #endif 14598fa52d88SHong Zhang ierr = PetscFree(rmap);CHKERRQ(ierr); 14608fa52d88SHong Zhang 14618fa52d88SHong Zhang for (i=0; i<ismax; i++) { 14628fa52d88SHong Zhang if (!allcolumns[i]) { 14638fa52d88SHong Zhang #if defined(PETSC_USE_CTABLE) 14648fa52d88SHong Zhang ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr); 14658fa52d88SHong Zhang #else 14668fa52d88SHong Zhang ierr = PetscFree(cmap[i]);CHKERRQ(ierr); 14678fa52d88SHong Zhang #endif 14688fa52d88SHong Zhang } 14698fa52d88SHong Zhang } 14708fa52d88SHong Zhang ierr = PetscFree(cmap);CHKERRQ(ierr); 1471606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1472d5b485abSSatish Balay 1473d5b485abSSatish Balay for (i=0; i<ismax; i++) { 147436f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147536f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1476d5b485abSSatish Balay } 14777a868f3eSHong Zhang 14787a868f3eSHong Zhang c->ijonly = PETSC_FALSE; /* set back to the default */ 14793a40ed3dSBarry Smith PetscFunctionReturn(0); 1480d5b485abSSatish Balay } 1481ca161407SBarry Smith 1482