123bdbc58SBarry Smith 2d5b485abSSatish Balay /* 3d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 4d5b485abSSatish Balay and to find submatrices that were shared across processors. 5d5b485abSSatish Balay */ 6c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h> 7c6db04a5SJed Brown #include <petscbt.h> 8d5b485abSSatish Balay 9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char**,PetscInt*,PetscInt**); 10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt**,PetscInt**,PetscInt*); 1109573ac7SBarry Smith extern PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 1209573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 13d5b485abSSatish Balay 144a2ae208SSatish Balay #undef __FUNCT__ 154a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ" 16b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov) 17d5b485abSSatish Balay { 186849ba73SBarry Smith PetscErrorCode ierr; 19d0f46423SBarry Smith PetscInt i,N=C->cmap->N, bs=C->rmap->bs; 2036f4e84dSSatish Balay IS *is_new; 2136f4e84dSSatish Balay 223a40ed3dSBarry Smith PetscFunctionBegin; 2382502324SSatish Balay ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr); 24df36731bSBarry Smith /* Convert the indices into block format */ 2505d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr); 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; 79*ce94432eSBarry 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 8724c049a4SHong Zhang ierr = PetscMalloc2(imax+1,const PetscInt*,&idx,imax,PetscInt,&n);CHKERRQ(ierr); 8824c049a4SHong Zhang 89d5b485abSSatish Balay for (i=0; i<imax; i++) { 90d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 91b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 92d5b485abSSatish Balay } 93d5b485abSSatish Balay 94d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 95d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 9623bdbc58SBarry Smith ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr); 9723bdbc58SBarry Smith ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 9823bdbc58SBarry Smith ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 9923bdbc58SBarry Smith ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr); 100d5b485abSSatish Balay for (i=0; i<imax; i++) { 101b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 102d5b485abSSatish Balay idx_i = idx[i]; 103d5b485abSSatish Balay len = n[i]; 104d5b485abSSatish Balay for (j=0; j<len; j++) { 105d5b485abSSatish Balay row = idx_i[j]; 106f23aa3ddSBarry Smith if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 107ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 108d5b485abSSatish Balay w4[proc]++; 109d5b485abSSatish Balay } 110d5b485abSSatish Balay for (j=0; j<size; j++) { 111d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 112d5b485abSSatish Balay } 113d5b485abSSatish Balay } 114d5b485abSSatish Balay 115d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 116d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1170e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 118d5b485abSSatish Balay w3[rank] = 0; 119d5b485abSSatish Balay for (i=0; i<size; i++) { 120d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 121d5b485abSSatish Balay } 122d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 123b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); 124d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 125d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 126d5b485abSSatish Balay } 127d5b485abSSatish Balay 128d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 129d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 130d5b485abSSatish Balay j = pa[i]; 131d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 132d5b485abSSatish Balay msz += w1[j]; 133d5b485abSSatish Balay } 134d5b485abSSatish Balay 135c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 136a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 137c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 138d5b485abSSatish Balay 139c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 140c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 141d5b485abSSatish Balay 142d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 14323bdbc58SBarry Smith ierr = PetscMalloc4(size,PetscInt*,&outdat,size,PetscInt*,&ptr,msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr); 14423bdbc58SBarry Smith ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr); 14523bdbc58SBarry Smith ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 146d5b485abSSatish Balay { 147b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 148d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 149d5b485abSSatish Balay j = pa[i]; 150d5b485abSSatish Balay iptr += ict; 151d5b485abSSatish Balay outdat[j] = iptr; 152d5b485abSSatish Balay ict = w1[j]; 153d5b485abSSatish Balay } 154d5b485abSSatish Balay } 155d5b485abSSatish Balay 156d5b485abSSatish Balay /* Form the outgoing messages */ 157d5b485abSSatish Balay /*plug in the headers*/ 158d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 159d5b485abSSatish Balay j = pa[i]; 160d5b485abSSatish Balay outdat[j][0] = 0; 161b24ad042SBarry Smith ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 162d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 163d5b485abSSatish Balay } 164d5b485abSSatish Balay 165d5b485abSSatish Balay /* Memory for doing local proc's work*/ 166d5b485abSSatish Balay { 1678f9f447aSHong Zhang ierr = PetscMalloc5(imax,PetscBT,&table, imax,PetscInt*,&data, imax,PetscInt,&isz, 1688f9f447aSHong Zhang Mbs*imax,PetscInt,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,char,&t_p);CHKERRQ(ierr); 1698f9f447aSHong Zhang ierr = PetscMemzero(table,imax*sizeof(PetscBT));CHKERRQ(ierr); 1708f9f447aSHong Zhang ierr = PetscMemzero(data,imax*sizeof(PetscInt*));CHKERRQ(ierr); 1718f9f447aSHong Zhang ierr = PetscMemzero(isz,imax*sizeof(PetscInt));CHKERRQ(ierr); 1728f9f447aSHong Zhang ierr = PetscMemzero(d_p,Mbs*imax*sizeof(PetscInt));CHKERRQ(ierr); 1738f9f447aSHong Zhang ierr = PetscMemzero(t_p,(Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char));CHKERRQ(ierr); 174d5b485abSSatish Balay 175bbd702dbSSatish Balay for (i=0; i<imax; i++) { 176b6410449SSatish Balay table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 177bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 178d5b485abSSatish Balay } 179d5b485abSSatish Balay } 180d5b485abSSatish Balay 181d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 182d5b485abSSatish Balay { 183b24ad042SBarry Smith PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 184f1af5d2fSBarry Smith PetscBT table_i; 185d5b485abSSatish Balay 186d5b485abSSatish Balay for (i=0; i<imax; i++) { 187b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 188d5b485abSSatish Balay n_i = n[i]; 189d5b485abSSatish Balay table_i = table[i]; 190d5b485abSSatish Balay idx_i = idx[i]; 191d5b485abSSatish Balay data_i = data[i]; 192d5b485abSSatish Balay isz_i = isz[i]; 193d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 194d5b485abSSatish Balay row = idx_i[j]; 195ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 196d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 197d5b485abSSatish Balay ctr[proc]++; 198d5b485abSSatish Balay *ptr[proc] = row; 199d5b485abSSatish Balay ptr[proc]++; 200d6b45a43SBarry Smith } else { /* Update the local table */ 20126fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 202d5b485abSSatish Balay } 203d5b485abSSatish Balay } 204d5b485abSSatish Balay /* Update the headers for the current IS */ 205d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 206d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 207d5b485abSSatish Balay outdat_j = outdat[j]; 208d5b485abSSatish Balay k = ++outdat_j[0]; 209d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 210d5b485abSSatish Balay outdat_j[2*k-1] = i; 211d5b485abSSatish Balay } 212d5b485abSSatish Balay } 213d5b485abSSatish Balay isz[i] = isz_i; 214d5b485abSSatish Balay } 215d5b485abSSatish Balay } 216d5b485abSSatish Balay 217d5b485abSSatish Balay /* Now post the sends */ 218b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 219d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 220d5b485abSSatish Balay j = pa[i]; 221b24ad042SBarry Smith ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 222d5b485abSSatish Balay } 223d5b485abSSatish Balay 224d5b485abSSatish Balay /* No longer need the original indices*/ 225d5b485abSSatish Balay for (i=0; i<imax; ++i) { 226d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 227d5b485abSSatish Balay } 22824c049a4SHong Zhang ierr = PetscFree2(idx,n);CHKERRQ(ierr); 229d5b485abSSatish Balay 230d5b485abSSatish Balay for (i=0; i<imax; ++i) { 2316bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 232d5b485abSSatish Balay } 233d5b485abSSatish Balay 234d5b485abSSatish Balay /* Do Local work*/ 235df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 236d5b485abSSatish Balay 237d5b485abSSatish Balay /* Receive messages*/ 238b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr); 2390c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);} 240d5b485abSSatish Balay 241b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 2420c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 243d5b485abSSatish Balay 244d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 24523bdbc58SBarry Smith ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr); 24623bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 247d5b485abSSatish Balay 248b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&xdata);CHKERRQ(ierr); 249b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt),&isz1);CHKERRQ(ierr); 250df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 251c05d87d6SBarry Smith ierr = PetscFree(rbuf[0]);CHKERRQ(ierr); 252606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 253d5b485abSSatish Balay 254d5b485abSSatish Balay /* Send the data back*/ 255d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 256d5b485abSSatish Balay { 257b24ad042SBarry Smith PetscMPIInt *rw1; 258d5b485abSSatish Balay 259b24ad042SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&rw1);CHKERRQ(ierr); 260b24ad042SBarry Smith ierr = PetscMemzero(rw1,size*sizeof(PetscInt));CHKERRQ(ierr); 261c7dd2924SSatish Balay 262d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 263d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 264e32f2f54SBarry Smith if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 265d5b485abSSatish Balay rw1[proc] = isz1[i]; 266d5b485abSSatish Balay } 267d5b485abSSatish Balay 268c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 269c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 270d5b485abSSatish Balay 271c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 272c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 27303512dc5SSatish Balay ierr = PetscFree(rw1);CHKERRQ(ierr); 274c7dd2924SSatish Balay } 275c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 276c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 277d5b485abSSatish Balay 278d5b485abSSatish Balay /* Now post the sends */ 279b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 280d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 281d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 282b24ad042SBarry Smith ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 283d5b485abSSatish Balay } 284d5b485abSSatish Balay 285d5b485abSSatish Balay /* receive work done on other processors*/ 286d5b485abSSatish Balay { 287b24ad042SBarry Smith PetscMPIInt idex; 288b24ad042SBarry Smith PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 289f1af5d2fSBarry Smith PetscBT table_i; 290d5b485abSSatish Balay MPI_Status *status2; 291d5b485abSSatish Balay 292169449a0SSatish Balay ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr); 293d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 294999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 295d5b485abSSatish Balay /* Process the message*/ 296999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 297d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 298999d9058SBarry Smith jmax = rbuf2[idex][0]; 299d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 300d5b485abSSatish Balay max = rbuf2_i[2*j]; 301d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 302d5b485abSSatish Balay isz_i = isz[is_no]; 303d5b485abSSatish Balay data_i = data[is_no]; 304d5b485abSSatish Balay table_i = table[is_no]; 305d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 306d5b485abSSatish Balay row = rbuf2_i[ct1]; 30726fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 308d5b485abSSatish Balay } 309d5b485abSSatish Balay isz[is_no] = isz_i; 310d5b485abSSatish Balay } 311d5b485abSSatish Balay } 3120c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);} 313606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 314d5b485abSSatish Balay } 315d5b485abSSatish Balay 316d5b485abSSatish Balay for (i=0; i<imax; ++i) { 31770b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 318d5b485abSSatish Balay } 319d5b485abSSatish Balay 320c7dd2924SSatish Balay 321c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 322c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 323c7dd2924SSatish Balay 324606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 325c05d87d6SBarry Smith ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr); 326606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 327606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 328606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 329606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 330606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 3318f9f447aSHong Zhang ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr); 332606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 333606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 334606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 335606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 336606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 3373a40ed3dSBarry Smith PetscFunctionReturn(0); 338d5b485abSSatish Balay } 339d5b485abSSatish Balay 3404a2ae208SSatish Balay #undef __FUNCT__ 3414a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local" 342d5b485abSSatish Balay /* 343df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 344d5b485abSSatish Balay the work on the local processor. 345d5b485abSSatish Balay 346d5b485abSSatish Balay Inputs: 347df36731bSBarry Smith C - MAT_MPIBAIJ; 348d5b485abSSatish Balay imax - total no of index sets processed at a time; 349df36731bSBarry Smith table - an array of char - size = Mbs bits. 350d5b485abSSatish Balay 351d5b485abSSatish Balay Output: 3520e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 353d5b485abSSatish Balay to each index set; 354d5b485abSSatish Balay data - pointer to the solutions 355d5b485abSSatish Balay */ 356b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 357d5b485abSSatish Balay { 358df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 359d5b485abSSatish Balay Mat A = c->A,B = c->B; 360df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 361b24ad042SBarry Smith PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 362b24ad042SBarry Smith PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 363f1af5d2fSBarry Smith PetscBT table_i; 364d5b485abSSatish Balay 3653a40ed3dSBarry Smith PetscFunctionBegin; 366899cda47SBarry Smith rstart = c->rstartbs; 367899cda47SBarry Smith cstart = c->cstartbs; 368d5b485abSSatish Balay ai = a->i; 369df36731bSBarry Smith aj = a->j; 370d5b485abSSatish Balay bi = b->i; 371df36731bSBarry Smith bj = b->j; 372d5b485abSSatish Balay garray = c->garray; 373d5b485abSSatish Balay 374d5b485abSSatish Balay 375d5b485abSSatish Balay for (i=0; i<imax; i++) { 376d5b485abSSatish Balay data_i = data[i]; 377d5b485abSSatish Balay table_i = table[i]; 378d5b485abSSatish Balay isz_i = isz[i]; 379d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 380d5b485abSSatish Balay row = data_i[j] - rstart; 381d5b485abSSatish Balay start = ai[row]; 382d5b485abSSatish Balay end = ai[row+1]; 383d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 384df36731bSBarry Smith val = aj[k] + cstart; 38526fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 386d5b485abSSatish Balay } 387d5b485abSSatish Balay start = bi[row]; 388d5b485abSSatish Balay end = bi[row+1]; 389d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 390df36731bSBarry Smith val = garray[bj[k]]; 39126fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 392d5b485abSSatish Balay } 393d5b485abSSatish Balay } 394d5b485abSSatish Balay isz[i] = isz_i; 395d5b485abSSatish Balay } 3963a40ed3dSBarry Smith PetscFunctionReturn(0); 397d5b485abSSatish Balay } 3984a2ae208SSatish Balay #undef __FUNCT__ 3994a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive" 400d5b485abSSatish Balay /* 401df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 402d5b485abSSatish Balay and return the output 403d5b485abSSatish Balay 404d5b485abSSatish Balay Input: 405d5b485abSSatish Balay C - the matrix 406d5b485abSSatish Balay nrqr - no of messages being processed. 407d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 408d5b485abSSatish Balay 409d5b485abSSatish Balay Output: 410d5b485abSSatish Balay xdata - array of messages to be sent back 411d5b485abSSatish Balay isz1 - size of each message 412d5b485abSSatish Balay 413a8c7a070SBarry Smith For better efficiency perhaps we should malloc separately each xdata[i], 414d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 4150e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of 416d5b485abSSatish Balay memory is used. 417d5b485abSSatish Balay 418d5b485abSSatish Balay */ 419b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 420d5b485abSSatish Balay { 421df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 422d5b485abSSatish Balay Mat A = c->A,B = c->B; 423df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 4246849ba73SBarry Smith PetscErrorCode ierr; 425b24ad042SBarry Smith PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 426b24ad042SBarry Smith PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 427d2a212eaSBarry Smith PetscInt val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 428b24ad042SBarry Smith PetscInt *rbuf_i,kmax,rbuf_0; 429f1af5d2fSBarry Smith PetscBT xtable; 430d5b485abSSatish Balay 4313a40ed3dSBarry Smith PetscFunctionBegin; 432df36731bSBarry Smith Mbs = c->Mbs; 433899cda47SBarry Smith rstart = c->rstartbs; 434899cda47SBarry Smith cstart = c->cstartbs; 435d5b485abSSatish Balay ai = a->i; 436df36731bSBarry Smith aj = a->j; 437d5b485abSSatish Balay bi = b->i; 438df36731bSBarry Smith bj = b->j; 439d5b485abSSatish Balay garray = c->garray; 440d5b485abSSatish Balay 441d5b485abSSatish Balay 442d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 443d5b485abSSatish Balay rbuf_i = rbuf[i]; 444d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 445d5b485abSSatish Balay ct += rbuf_0; 44626fbe8dcSKarl Rupp for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j]; 447d5b485abSSatish Balay } 448d5b485abSSatish Balay 449701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 450701b8814SBarry Smith else max1 = 1; 451d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 452b24ad042SBarry Smith ierr = PetscMalloc(mem_estimate*sizeof(PetscInt),&xdata[0]);CHKERRQ(ierr); 453d5b485abSSatish Balay ++no_malloc; 45453b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr); 455b24ad042SBarry Smith ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr); 456d5b485abSSatish Balay 457d5b485abSSatish Balay ct3 = 0; 458d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 459d5b485abSSatish Balay rbuf_i = rbuf[i]; 460d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 461d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 462d5b485abSSatish Balay ct2 = ct1; 463d5b485abSSatish Balay ct3 += ct1; 464d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 4656831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 466d5b485abSSatish Balay oct2 = ct2; 467d5b485abSSatish Balay kmax = rbuf_i[2*j]; 468d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 469d5b485abSSatish Balay row = rbuf_i[ct1]; 470f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 471d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 472b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 473b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 474b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 475606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 476d5b485abSSatish Balay xdata[0] = tmp; 477d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 47826fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 479d5b485abSSatish Balay } 480d5b485abSSatish Balay xdata[i][ct2++] = row; 481d5b485abSSatish Balay ct3++; 482d5b485abSSatish Balay } 483d5b485abSSatish Balay } 484d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 485d5b485abSSatish Balay row = xdata[i][k] - rstart; 486d5b485abSSatish Balay start = ai[row]; 487d5b485abSSatish Balay end = ai[row+1]; 488d5b485abSSatish Balay for (l=start; l<end; l++) { 489df36731bSBarry Smith val = aj[l] + cstart; 490f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 491d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 492b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 493b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 494b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 495606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 496d5b485abSSatish Balay xdata[0] = tmp; 497d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 49826fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 499d5b485abSSatish Balay } 500d5b485abSSatish Balay xdata[i][ct2++] = val; 501d5b485abSSatish Balay ct3++; 502d5b485abSSatish Balay } 503d5b485abSSatish Balay } 504d5b485abSSatish Balay start = bi[row]; 505d5b485abSSatish Balay end = bi[row+1]; 506d5b485abSSatish Balay for (l=start; l<end; l++) { 507df36731bSBarry Smith val = garray[bj[l]]; 508f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 509d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 510b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 511b24ad042SBarry Smith ierr = PetscMalloc(new_estimate * sizeof(PetscInt),&tmp);CHKERRQ(ierr); 512b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 513606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 514d5b485abSSatish Balay xdata[0] = tmp; 515d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 51626fbe8dcSKarl Rupp for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 517d5b485abSSatish Balay } 518d5b485abSSatish Balay xdata[i][ct2++] = val; 519d5b485abSSatish Balay ct3++; 520d5b485abSSatish Balay } 521d5b485abSSatish Balay } 522d5b485abSSatish Balay } 523d5b485abSSatish Balay /* Update the header*/ 524d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 525d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 526d5b485abSSatish Balay } 527d5b485abSSatish Balay xdata[i][0] = rbuf_0; 528d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 529d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 530d5b485abSSatish Balay } 53194bacf5dSBarry Smith ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr); 5321e2582c4SBarry Smith ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr); 5333a40ed3dSBarry Smith PetscFunctionReturn(0); 534d5b485abSSatish Balay } 535d5b485abSSatish Balay 5364a2ae208SSatish Balay #undef __FUNCT__ 5374a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ" 538b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 539d5b485abSSatish Balay { 54036f4e84dSSatish Balay IS *isrow_new,*iscol_new; 541cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 5426849ba73SBarry Smith PetscErrorCode ierr; 5434da72fa9SHong Zhang PetscInt nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs; 5444da72fa9SHong Zhang PetscBool colflag,*allcolumns,*allrows; 545a2feddc5SSatish Balay 5463a40ed3dSBarry Smith PetscFunctionBegin; 54729dcf524SDmitry Karpeev /* Currently, unsorted column indices will result in inverted column indices in the resulting submatrices. */ 54829dcf524SDmitry Karpeev for (i = 0; i < ismax; ++i) { 54929dcf524SDmitry Karpeev PetscBool sorted; 55029dcf524SDmitry Karpeev ierr = ISSorted(iscol[i], &sorted);CHKERRQ(ierr); 551c7e6e2c7SJed Brown if (!sorted) SETERRQ1(((PetscObject)iscol[i])->comm, PETSC_ERR_SUP, "Column index set %D not sorted", i); 55229dcf524SDmitry Karpeev } 55329dcf524SDmitry Karpeev /* The compression and expansion should be avoided. Doesn't point 55429dcf524SDmitry Karpeev out errors, might change the indices, hence buggey */ 55523bdbc58SBarry Smith ierr = PetscMalloc2(ismax+1,IS,&isrow_new,ismax+1,IS,&iscol_new);CHKERRQ(ierr); 55605d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_new);CHKERRQ(ierr); 55705d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_new);CHKERRQ(ierr); 558cf4f527aSSatish Balay 559307b7a18SHong Zhang /* Check for special case: each processor gets entire matrix columns */ 5604da72fa9SHong Zhang ierr = PetscMalloc2(ismax+1,PetscBool,&allcolumns,ismax+1,PetscBool,&allrows);CHKERRQ(ierr); 561307b7a18SHong Zhang for (i=0; i<ismax; i++) { 562307b7a18SHong Zhang ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 563307b7a18SHong Zhang ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr); 56426fbe8dcSKarl Rupp if (colflag && ncol == C->cmap->N) allcolumns[i] = PETSC_TRUE; 56526fbe8dcSKarl Rupp else allcolumns[i] = PETSC_FALSE; 5664da72fa9SHong Zhang 5674da72fa9SHong Zhang ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 5684da72fa9SHong Zhang ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr); 56926fbe8dcSKarl Rupp if (colflag && nrow == C->rmap->N) allrows[i] = PETSC_TRUE; 57026fbe8dcSKarl Rupp else allrows[i] = PETSC_FALSE; 571307b7a18SHong Zhang } 572307b7a18SHong Zhang 573cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 574cf4f527aSSatish Balay if (scall != MAT_REUSE_MATRIX) { 57582502324SSatish Balay ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr); 576cf4f527aSSatish Balay } 577cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 578b24ad042SBarry Smith nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 579329f5518SBarry Smith if (!nmax) nmax = 1; 580cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); 581cf4f527aSSatish Balay 582653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 583*ce94432eSBarry Smith ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 584cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 585cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 586cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 587cf4f527aSSatish Balay else max_no = ismax-pos; 5884da72fa9SHong Zhang ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 589cf4f527aSSatish Balay pos += max_no; 590cf4f527aSSatish Balay } 59136f4e84dSSatish Balay 59236f4e84dSSatish Balay for (i=0; i<ismax; i++) { 5936bf464f9SBarry Smith ierr = ISDestroy(&isrow_new[i]);CHKERRQ(ierr); 5946bf464f9SBarry Smith ierr = ISDestroy(&iscol_new[i]);CHKERRQ(ierr); 59536f4e84dSSatish Balay } 59623bdbc58SBarry Smith ierr = PetscFree2(isrow_new,iscol_new);CHKERRQ(ierr); 5974da72fa9SHong Zhang ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr); 5983a40ed3dSBarry Smith PetscFunctionReturn(0); 599a2feddc5SSatish Balay } 600a2feddc5SSatish Balay 601233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 6024a2ae208SSatish Balay #undef __FUNCT__ 6034a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc" 604e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 605233c2ff6SSatish Balay { 606e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 6074dc2109aSBarry Smith PetscMPIInt fproc; 6084dc2109aSBarry Smith PetscErrorCode ierr; 609233c2ff6SSatish Balay 610233c2ff6SSatish Balay PetscFunctionBegin; 6114dc2109aSBarry Smith ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr); 61223ce1328SBarry Smith if (fproc > size) fproc = size; 613e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 614e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 615233c2ff6SSatish Balay else fproc++; 616233c2ff6SSatish Balay } 617e005ede5SBarry Smith *rank = fproc; 618233c2ff6SSatish Balay PetscFunctionReturn(0); 619233c2ff6SSatish Balay } 620233c2ff6SSatish Balay #endif 621233c2ff6SSatish Balay 622a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 623b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 6244a2ae208SSatish Balay #undef __FUNCT__ 6254a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local" 6264da72fa9SHong Zhang PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 627a2feddc5SSatish Balay { 628df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 629cf4f527aSSatish Balay Mat A = c->A; 630df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 6315d0c19d7SBarry Smith const PetscInt **irow,**icol,*irow_i; 6325d0c19d7SBarry Smith PetscInt *nrow,*ncol,*w3,*w4,start; 6336849ba73SBarry Smith PetscErrorCode ierr; 6349405f653SBarry Smith PetscMPIInt size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc; 6359405f653SBarry Smith PetscInt **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row; 636b24ad042SBarry Smith PetscInt nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol; 637b24ad042SBarry Smith PetscInt **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 6385d0c19d7SBarry Smith PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 639052f0c41SBarry Smith PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 640d0f46423SBarry Smith PetscInt bs =C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 641899cda47SBarry Smith PetscInt cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 642899cda47SBarry Smith PetscInt *bmap = c->garray,ctmp,rstart=c->rstartbs; 6437a868f3eSHong Zhang MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3; 6447a868f3eSHong Zhang MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3; 645d5b485abSSatish Balay MPI_Comm comm; 646ace3abfcSBarry Smith PetscBool flag; 647b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1; 6487a868f3eSHong Zhang PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 64926fbe8dcSKarl Rupp 6507a868f3eSHong Zhang /* variables below are used for the matrix numerical values - case of !ijonly */ 6517a868f3eSHong Zhang MPI_Request *r_waits4,*s_waits4; 6527a868f3eSHong Zhang MPI_Status *r_status4,*s_status4; 6530298fd71SBarry Smith MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a = NULL,*sbuf_aa_i,*vworkA = NULL,*vworkB = NULL; 6547a868f3eSHong Zhang MatScalar *a_a=a->a,*b_a=b->a; 655c7dd2924SSatish Balay 65680d608c0SSatish Balay #if defined(PETSC_USE_CTABLE) 657b24ad042SBarry Smith PetscInt tt; 6580298fd71SBarry Smith PetscTable *rmap,*cmap,rmap_i,cmap_i=NULL; 659233c2ff6SSatish Balay #else 6600298fd71SBarry Smith PetscInt **cmap,*cmap_i=NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 661233c2ff6SSatish Balay #endif 662d5b485abSSatish Balay 6633a40ed3dSBarry Smith PetscFunctionBegin; 664*ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 6657adad957SLisandro Dalcin tag0 = ((PetscObject)C)->tag; 666d5b485abSSatish Balay size = c->size; 667d5b485abSSatish Balay rank = c->rank; 668d5b485abSSatish Balay 66934e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 67034e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 67134e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 67234e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 67334e6ae18SSatish Balay 674052f0c41SBarry Smith #if defined(PETSC_USE_CTABLE) 67523bdbc58SBarry Smith ierr = PetscMalloc4(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol);CHKERRQ(ierr); 676052f0c41SBarry Smith #else 67723bdbc58SBarry Smith ierr = PetscMalloc5(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol,Mbs+1,PetscInt,&rtable);CHKERRQ(ierr); 678d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 679d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 68082a7e548SBarry Smith jmax = C->rmap->range[i+1]/bs; 68126fbe8dcSKarl Rupp for (; j<jmax; j++) rtable[j] = i; 682d5b485abSSatish Balay } 683233c2ff6SSatish Balay #endif 684233c2ff6SSatish Balay 685233c2ff6SSatish Balay for (i=0; i<ismax; i++) { 6864da72fa9SHong Zhang if (allrows[i]) { 6870298fd71SBarry Smith irow[i] = NULL; 6884da72fa9SHong Zhang nrow[i] = C->rmap->N/bs; 6894da72fa9SHong Zhang } else { 690233c2ff6SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 691233c2ff6SSatish Balay ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 6924da72fa9SHong Zhang } 6934da72fa9SHong Zhang 694307b7a18SHong Zhang if (allcolumns[i]) { 6950298fd71SBarry Smith icol[i] = NULL; 696307b7a18SHong Zhang ncol[i] = C->cmap->N/bs; 697307b7a18SHong Zhang } else { 698307b7a18SHong Zhang ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 699233c2ff6SSatish Balay ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 700233c2ff6SSatish Balay } 701307b7a18SHong Zhang } 702d5b485abSSatish Balay 703d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg,and buffer space 704d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 70523bdbc58SBarry Smith ierr = PetscMalloc4(size,PetscMPIInt,&w1,size,PetscMPIInt,&w2,size,PetscInt,&w3,size,PetscInt,&w4);CHKERRQ(ierr); 70623bdbc58SBarry Smith ierr = PetscMemzero(w1,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 70723bdbc58SBarry Smith ierr = PetscMemzero(w2,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 70823bdbc58SBarry Smith ierr = PetscMemzero(w3,size*sizeof(PetscInt));CHKERRQ(ierr); 709d5b485abSSatish Balay for (i=0; i<ismax; i++) { 710b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 711d5b485abSSatish Balay jmax = nrow[i]; 712d5b485abSSatish Balay irow_i = irow[i]; 713d5b485abSSatish Balay for (j=0; j<jmax; j++) { 71426fbe8dcSKarl Rupp if (allrows[i]) row = j; 71526fbe8dcSKarl Rupp else row = irow_i[j]; 71626fbe8dcSKarl Rupp 717233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 718899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&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 */ 75923bdbc58SBarry Smith ierr = PetscMalloc4(size,PetscInt*,&sbuf1,size,PetscInt*,&ptr,2*msz,PetscInt,&tmp,size,PetscInt,&ctr);CHKERRQ(ierr); 76023bdbc58SBarry Smith ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 76123bdbc58SBarry Smith ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 762d5b485abSSatish Balay { 763b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 764d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 765d5b485abSSatish Balay j = pa[i]; 766d5b485abSSatish Balay iptr += ict; 767d5b485abSSatish Balay sbuf1[j] = iptr; 768d5b485abSSatish Balay ict = w1[j]; 769d5b485abSSatish Balay } 770d5b485abSSatish Balay } 771d5b485abSSatish Balay 772d5b485abSSatish Balay /* Form the outgoing messages */ 773d5b485abSSatish Balay /* Initialise the header space */ 774d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 775d5b485abSSatish Balay j = pa[i]; 776d5b485abSSatish Balay sbuf1[j][0] = 0; 777b24ad042SBarry Smith ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 778d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 779d5b485abSSatish Balay } 780d5b485abSSatish Balay 781d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 782d5b485abSSatish Balay for (i=0; i<ismax; i++) { 783b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 784d5b485abSSatish Balay irow_i = irow[i]; 785d5b485abSSatish Balay jmax = nrow[i]; 786d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 78726fbe8dcSKarl Rupp if (allrows[i]) row = j; 78826fbe8dcSKarl Rupp else row = irow_i[j]; 78926fbe8dcSKarl Rupp 790233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 791899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 792233c2ff6SSatish Balay #else 793d5b485abSSatish Balay proc = rtable[row]; 794233c2ff6SSatish Balay #endif 795d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 796d5b485abSSatish Balay ctr[proc]++; 797d5b485abSSatish Balay *ptr[proc] = row; 798d5b485abSSatish Balay ptr[proc]++; 799d5b485abSSatish Balay } 800d5b485abSSatish Balay } 801d5b485abSSatish Balay /* Update the headers for the current IS */ 802d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 80306ef35abSSatish Balay if ((ctr_j = ctr[j])) { 804d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 805d5b485abSSatish Balay k = ++sbuf1_j[0]; 806d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 807d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 808d5b485abSSatish Balay } 809d5b485abSSatish Balay } 810d5b485abSSatish Balay } 811d5b485abSSatish Balay 812d5b485abSSatish Balay /* Now post the sends */ 813b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 814d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 815d5b485abSSatish Balay j = pa[i]; 816b24ad042SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 817d5b485abSSatish Balay } 818d5b485abSSatish Balay 819d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 820b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr); 821b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf2);CHKERRQ(ierr); 822d5b485abSSatish Balay rbuf2[0] = tmp + msz; 823d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 824d5b485abSSatish Balay j = pa[i]; 825d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 826d5b485abSSatish Balay } 827d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 828d5b485abSSatish Balay j = pa[i]; 829b24ad042SBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 830d5b485abSSatish Balay } 831d5b485abSSatish Balay 832d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 833d5b485abSSatish Balay 834d5b485abSSatish Balay /* Receive messages*/ 835b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 836b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr); 83723bdbc58SBarry Smith ierr = PetscMalloc3(nrqr+1,PetscInt*,&sbuf2,nrqr,PetscInt,&req_size,nrqr,PetscInt,&req_source);CHKERRQ(ierr); 838d5b485abSSatish Balay { 839df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 840b24ad042SBarry Smith PetscInt *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 841d5b485abSSatish Balay 842d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 843999d9058SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 84426fbe8dcSKarl Rupp 845999d9058SBarry Smith req_size[idex] = 0; 846999d9058SBarry Smith rbuf1_i = rbuf1[idex]; 847d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 848b24ad042SBarry Smith ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 849b24ad042SBarry Smith ierr = PetscMalloc(end*sizeof(PetscInt),&sbuf2[idex]);CHKERRQ(ierr); 850999d9058SBarry Smith sbuf2_i = sbuf2[idex]; 851d5b485abSSatish Balay for (j=start; j<end; j++) { 852d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 853d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 854d5b485abSSatish Balay sbuf2_i[j] = ncols; 855999d9058SBarry Smith req_size[idex] += ncols; 856d5b485abSSatish Balay } 857999d9058SBarry Smith req_source[idex] = r_status1[i].MPI_SOURCE; 858d5b485abSSatish Balay /* form the header */ 859999d9058SBarry Smith sbuf2_i[0] = req_size[idex]; 86026fbe8dcSKarl Rupp for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 861b24ad042SBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 862d5b485abSSatish Balay } 863d5b485abSSatish Balay } 864606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 865606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 866d5b485abSSatish Balay 867d5b485abSSatish Balay /* recv buffer sizes */ 868d5b485abSSatish Balay /* Receive messages*/ 869b24ad042SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt*),&rbuf3);CHKERRQ(ierr); 870b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr); 871b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr); 8727a868f3eSHong Zhang if (!ijonly) { 8737a868f3eSHong Zhang ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr); 8747a868f3eSHong Zhang ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr); 8757a868f3eSHong Zhang } 876d5b485abSSatish Balay 877d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 878999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 879b24ad042SBarry Smith ierr = PetscMalloc(rbuf2[idex][0]*sizeof(PetscInt),&rbuf3[idex]);CHKERRQ(ierr); 880b24ad042SBarry Smith ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 8817a868f3eSHong Zhang if (!ijonly) { 8827a868f3eSHong Zhang ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr); 883b24ad042SBarry Smith ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 884d5b485abSSatish Balay } 8857a868f3eSHong Zhang } 886606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 887606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 888d5b485abSSatish Balay 889d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 890b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr); 891b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr); 892d5b485abSSatish Balay 8930c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 8940c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 895606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 896606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 897606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 898606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 899d5b485abSSatish Balay 900d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 901b24ad042SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscInt*),&sbuf_aj);CHKERRQ(ierr); 902d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 903b24ad042SBarry Smith ierr = PetscMalloc((j+1)*sizeof(PetscInt),&sbuf_aj[0]);CHKERRQ(ierr); 904d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 905d5b485abSSatish Balay 906b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr); 907d5b485abSSatish Balay { 908d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 909d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 910d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 911d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 912d5b485abSSatish Balay ct2 = 0; 913d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 914d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 915d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 916e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 917d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 918d5b485abSSatish Balay ncols = nzA + nzB; 919d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 920d5b485abSSatish Balay 921d5b485abSSatish Balay /* load the column indices for this row into cols*/ 922d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 923d5b485abSSatish Balay for (l=0; l<nzB; l++) { 924d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 925d5b485abSSatish Balay else break; 926d5b485abSSatish Balay } 927d5b485abSSatish Balay imark = l; 928d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 929d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 930d5b485abSSatish Balay ct2 += ncols; 931d5b485abSSatish Balay } 932d5b485abSSatish Balay } 933b24ad042SBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 934d5b485abSSatish Balay } 935d5b485abSSatish Balay } 936b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr); 937b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr); 938d5b485abSSatish Balay 939d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 9407a868f3eSHong Zhang if (!ijonly) { 94182502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar*),&sbuf_aa);CHKERRQ(ierr); 942d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 94382502324SSatish Balay ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr); 944a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 945d5b485abSSatish Balay 946b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr); 947d5b485abSSatish Balay { 948d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 949d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 950d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 951d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 952d5b485abSSatish Balay ct2 = 0; 953d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 954d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 955d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 956e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 957d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 958d5b485abSSatish Balay ncols = nzA + nzB; 959d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 960a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 961d5b485abSSatish Balay 962d5b485abSSatish Balay /* load the column values for this row into vals*/ 9635b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 964d5b485abSSatish Balay for (l=0; l<nzB; l++) { 9653a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 9663eda8832SBarry Smith ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 96726fbe8dcSKarl Rupp } else break; 968d5b485abSSatish Balay } 969d5b485abSSatish Balay imark = l; 9703a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 9713eda8832SBarry Smith ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9723a40ed3dSBarry Smith } 9733a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 9743eda8832SBarry Smith ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 9753a40ed3dSBarry Smith } 976d5b485abSSatish Balay ct2 += ncols; 977d5b485abSSatish Balay } 978d5b485abSSatish Balay } 9793eda8832SBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 980d5b485abSSatish Balay } 981d5b485abSSatish Balay } 982b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr); 983b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr); 9847a868f3eSHong Zhang } 985533163c2SBarry Smith ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 986606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 987d5b485abSSatish Balay 988d5b485abSSatish Balay /* Form the matrix */ 989307b7a18SHong Zhang /* create col map: global col of C -> local col of submatrices */ 990d5b485abSSatish Balay { 9915d0c19d7SBarry Smith const PetscInt *icol_i; 992233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 9938fa52d88SHong Zhang ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&cmap);CHKERRQ(ierr); 994ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 995307b7a18SHong Zhang if (!allcolumns[i]) { 9968fa52d88SHong Zhang ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 997307b7a18SHong Zhang jmax = ncol[i]; 998307b7a18SHong Zhang icol_i = icol[i]; 9998fa52d88SHong Zhang cmap_i = cmap[i]; 1000307b7a18SHong Zhang for (j=0; j<jmax; j++) { 10013861aac3SJed Brown ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1002307b7a18SHong Zhang } 1003307b7a18SHong Zhang } else { 10040298fd71SBarry Smith cmap[i] = NULL; 1005307b7a18SHong Zhang } 1006233c2ff6SSatish Balay } 1007233c2ff6SSatish Balay #else 100823bdbc58SBarry Smith ierr = PetscMalloc(ismax*sizeof(PetscInt*),&cmap);CHKERRQ(ierr); 1009d5b485abSSatish Balay for (i=0; i<ismax; i++) { 10108fa52d88SHong Zhang if (!allcolumns[i]) { 10118fa52d88SHong Zhang ierr = PetscMalloc(c->Nbs*sizeof(PetscInt),&cmap[i]);CHKERRQ(ierr); 10128fa52d88SHong Zhang ierr = PetscMemzero(cmap[i],c->Nbs*sizeof(PetscInt));CHKERRQ(ierr); 1013d5b485abSSatish Balay jmax = ncol[i]; 1014d5b485abSSatish Balay icol_i = icol[i]; 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 } 10198fa52d88SHong Zhang } else { /* allcolumns[i] */ 10200298fd71SBarry Smith cmap[i] = NULL; 10218fa52d88SHong Zhang } 1022d5b485abSSatish Balay } 1023307b7a18SHong Zhang #endif 1024d5b485abSSatish Balay } 1025d5b485abSSatish Balay 1026d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 102726fbe8dcSKarl Rupp for (i=0,j=0; i<ismax; i++) j += nrow[i]; 1028052f0c41SBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr); 1029b24ad042SBarry Smith lens[0] = (PetscInt*)(lens + ismax); 1030b24ad042SBarry Smith ierr = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr); 103126fbe8dcSKarl Rupp for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 1032d5b485abSSatish Balay 1033d5b485abSSatish Balay /* Update lens from local data */ 1034d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1035d5b485abSSatish Balay jmax = nrow[i]; 10368fa52d88SHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 1037d5b485abSSatish Balay irow_i = irow[i]; 1038d5b485abSSatish Balay lens_i = lens[i]; 1039d5b485abSSatish Balay for (j=0; j<jmax; j++) { 104026fbe8dcSKarl Rupp if (allrows[i]) row = j; 104126fbe8dcSKarl Rupp else row = irow_i[j]; 104226fbe8dcSKarl Rupp 1043233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1044899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&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]; 1053307b7a18SHong Zhang if (!allcolumns[i]) { 1054233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1055233c2ff6SSatish Balay for (k=0; k<nzA; k++) { 10568fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 105726fbe8dcSKarl Rupp if (tt) lens_i[j]++; 1058233c2ff6SSatish Balay } 1059233c2ff6SSatish Balay for (k=0; k<nzB; k++) { 10608fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 106126fbe8dcSKarl Rupp if (tt) lens_i[j]++; 1062233c2ff6SSatish Balay } 1063307b7a18SHong Zhang 1064233c2ff6SSatish Balay #else 1065ca161407SBarry Smith for (k=0; k<nzA; k++) { 106626fbe8dcSKarl Rupp if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1067ca161407SBarry Smith } 1068ca161407SBarry Smith for (k=0; k<nzB; k++) { 106926fbe8dcSKarl Rupp if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1070d5b485abSSatish Balay } 1071233c2ff6SSatish Balay #endif 1072307b7a18SHong Zhang } else { /* allcolumns */ 1073307b7a18SHong Zhang lens_i[j] = nzA + nzB; 1074307b7a18SHong Zhang } 1075d5b485abSSatish Balay } 1076a2feddc5SSatish Balay } 1077ca161407SBarry Smith } 1078233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1079d5b485abSSatish Balay /* Create row map*/ 10808fa52d88SHong Zhang ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rmap);CHKERRQ(ierr); 1081ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 10828fa52d88SHong Zhang ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1083233c2ff6SSatish Balay } 1084233c2ff6SSatish Balay #else 1085233c2ff6SSatish Balay /* Create row map*/ 1086052f0c41SBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr); 1087b24ad042SBarry Smith rmap[0] = (PetscInt*)(rmap + ismax); 1088b24ad042SBarry Smith ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr); 108926fbe8dcSKarl Rupp for (i=1; i<ismax; i++) rmap[i] = rmap[i-1] + Mbs; 1090233c2ff6SSatish Balay #endif 1091d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1092d5b485abSSatish Balay irow_i = irow[i]; 1093d5b485abSSatish Balay jmax = nrow[i]; 1094233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 10958fa52d88SHong Zhang rmap_i = rmap[i]; 1096233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 10974da72fa9SHong Zhang if (allrows[i]) { 10983861aac3SJed Brown ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 10994da72fa9SHong Zhang } else { 11003861aac3SJed Brown ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1101233c2ff6SSatish Balay } 11024da72fa9SHong Zhang } 1103233c2ff6SSatish Balay #else 1104233c2ff6SSatish Balay rmap_i = rmap[i]; 1105d5b485abSSatish Balay for (j=0; j<jmax; j++) { 110626fbe8dcSKarl Rupp if (allrows[i]) rmap_i[j] = j; 110726fbe8dcSKarl Rupp else rmap_i[irow_i[j]] = j; 11084da72fa9SHong Zhang } 1109233c2ff6SSatish Balay #endif 1110d5b485abSSatish Balay } 1111d5b485abSSatish Balay 1112d5b485abSSatish Balay /* Update lens from offproc data */ 1113d5b485abSSatish Balay { 1114b24ad042SBarry Smith PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1115b24ad042SBarry Smith PetscMPIInt ii; 1116d5b485abSSatish Balay 1117d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1118b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr); 1119b24ad042SBarry Smith idex = pa[ii]; 1120999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1121d5b485abSSatish Balay jmax = sbuf1_i[0]; 1122d5b485abSSatish Balay ct1 = 2*jmax+1; 1123d5b485abSSatish Balay ct2 = 0; 1124b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1125b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 1126d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1127d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1128d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1129d5b485abSSatish Balay lens_i = lens[is_no]; 11308fa52d88SHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1131d5b485abSSatish Balay rmap_i = rmap[is_no]; 1132d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1133233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 11348fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1135233c2ff6SSatish Balay row--; 113626fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1137233c2ff6SSatish Balay #else 1138d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1139233c2ff6SSatish Balay #endif 1140d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1141d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1142307b7a18SHong Zhang if (!allcolumns[is_no]) { 1143233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 11448fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 114526fbe8dcSKarl Rupp if (tt) lens_i[row]++; 1146233c2ff6SSatish Balay #else 114726fbe8dcSKarl Rupp if (cmap_i[rbuf3_i[ct2]]) lens_i[row]++; 1148233c2ff6SSatish Balay #endif 1149307b7a18SHong Zhang } else { /* allcolumns */ 1150307b7a18SHong Zhang lens_i[row]++; 1151307b7a18SHong Zhang } 1152d5b485abSSatish Balay } 1153d5b485abSSatish Balay } 1154d5b485abSSatish Balay } 1155d5b485abSSatish Balay } 1156d5b485abSSatish Balay } 1157606d414cSSatish Balay ierr = PetscFree(r_status3);CHKERRQ(ierr); 1158606d414cSSatish Balay ierr = PetscFree(r_waits3);CHKERRQ(ierr); 11590c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 1160606d414cSSatish Balay ierr = PetscFree(s_status3);CHKERRQ(ierr); 1161606d414cSSatish Balay ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1162d5b485abSSatish Balay 1163d5b485abSSatish Balay /* Create the submatrices */ 1164d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 11657a868f3eSHong Zhang if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet"); 1166d5b485abSSatish Balay /* 1167d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1168d5b485abSSatish Balay */ 1169d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1170df36731bSBarry Smith mat = (Mat_SeqBAIJ*)(submats[i]->data); 1171e7e72b3dSBarry 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"); 1172b24ad042SBarry Smith ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 1173e7e72b3dSBarry Smith if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 1174d5b485abSSatish Balay /* Initial matrix as if empty */ 1175b24ad042SBarry Smith ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr); 117626fbe8dcSKarl Rupp 1177d5f3da31SBarry Smith submats[i]->factortype = C->factortype; 1178d5b485abSSatish Balay } 1179ca161407SBarry Smith } else { 11807a868f3eSHong Zhang PetscInt bs_tmp; 118126fbe8dcSKarl Rupp if (ijonly) bs_tmp = 1; 118226fbe8dcSKarl Rupp else bs_tmp = bs; 1183d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1184f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 11857a868f3eSHong Zhang ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr); 11867adad957SLisandro Dalcin ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 11877a868f3eSHong Zhang ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 11887a868f3eSHong Zhang ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 1189d5b485abSSatish Balay } 1190d5b485abSSatish Balay } 1191d5b485abSSatish Balay 1192d5b485abSSatish Balay /* Assemble the matrices */ 1193d5b485abSSatish Balay /* First assemble the local rows */ 1194d5b485abSSatish Balay { 1195b24ad042SBarry Smith PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 11960298fd71SBarry Smith MatScalar *imat_a = NULL; 1197d5b485abSSatish Balay 1198d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1199df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[i]->data; 1200d5b485abSSatish Balay imat_ilen = mat->ilen; 1201d5b485abSSatish Balay imat_j = mat->j; 1202d5b485abSSatish Balay imat_i = mat->i; 12037a868f3eSHong Zhang if (!ijonly) imat_a = mat->a; 12048fa52d88SHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 1205d5b485abSSatish Balay rmap_i = rmap[i]; 1206d5b485abSSatish Balay irow_i = irow[i]; 1207d5b485abSSatish Balay jmax = nrow[i]; 1208d5b485abSSatish Balay for (j=0; j<jmax; j++) { 120926fbe8dcSKarl Rupp if (allrows[i]) row = j; 121026fbe8dcSKarl Rupp else row = irow_i[j]; 121126fbe8dcSKarl Rupp 1212233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1213899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&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]; 12237a868f3eSHong Zhang if (!ijonly) { 122436f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 122536f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 12267a868f3eSHong Zhang } 1227233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 12288fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1229233c2ff6SSatish Balay row--; 123026fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1231233c2ff6SSatish Balay #else 123236f4e84dSSatish Balay row = rmap_i[row + rstart]; 1233233c2ff6SSatish Balay #endif 1234df36731bSBarry Smith mat_i = imat_i[row]; 12357a868f3eSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1236d5b485abSSatish Balay mat_j = imat_j + mat_i; 123736f4e84dSSatish Balay ilen_row = imat_ilen[row]; 123836f4e84dSSatish Balay 123936f4e84dSSatish Balay /* load the column indices for this row into cols*/ 1240307b7a18SHong Zhang if (!allcolumns[i]) { 124136f4e84dSSatish Balay for (l=0; l<nzB; l++) { 124236f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 1243233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 12448fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1245233c2ff6SSatish Balay if (tcol) { 1246233c2ff6SSatish Balay #else 124736f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1248233c2ff6SSatish Balay #endif 1249df36731bSBarry Smith *mat_j++ = tcol - 1; 12503eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1251549d3d68SSatish Balay mat_a += bs2; 1252d5b485abSSatish Balay ilen_row++; 1253d5b485abSSatish Balay } 1254ca161407SBarry Smith } else break; 125536f4e84dSSatish Balay } 125636f4e84dSSatish Balay imark = l; 125736f4e84dSSatish Balay for (l=0; l<nzA; l++) { 1258233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 12598fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1260233c2ff6SSatish Balay if (tcol) { 1261233c2ff6SSatish Balay #else 126236f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 1263233c2ff6SSatish Balay #endif 126436f4e84dSSatish Balay *mat_j++ = tcol - 1; 12657a868f3eSHong Zhang if (!ijonly) { 12663eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1267549d3d68SSatish Balay mat_a += bs2; 12687a868f3eSHong Zhang } 126936f4e84dSSatish Balay ilen_row++; 127036f4e84dSSatish Balay } 127136f4e84dSSatish Balay } 127236f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 1273233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 12748fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1275233c2ff6SSatish Balay if (tcol) { 1276233c2ff6SSatish Balay #else 127736f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1278233c2ff6SSatish Balay #endif 127936f4e84dSSatish Balay *mat_j++ = tcol - 1; 12807a868f3eSHong Zhang if (!ijonly) { 12813eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1282549d3d68SSatish Balay mat_a += bs2; 12837a868f3eSHong Zhang } 128436f4e84dSSatish Balay ilen_row++; 128536f4e84dSSatish Balay } 128636f4e84dSSatish Balay } 1287307b7a18SHong Zhang } else { /* allcolumns */ 1288307b7a18SHong Zhang for (l=0; l<nzB; l++) { 1289307b7a18SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1290307b7a18SHong Zhang *mat_j++ = ctmp; 1291307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1292307b7a18SHong Zhang mat_a += bs2; 1293307b7a18SHong Zhang ilen_row++; 1294307b7a18SHong Zhang } else break; 1295307b7a18SHong Zhang } 1296307b7a18SHong Zhang imark = l; 1297307b7a18SHong Zhang for (l=0; l<nzA; l++) { 1298307b7a18SHong Zhang *mat_j++ = cstart+cworkA[l]; 1299307b7a18SHong Zhang if (!ijonly) { 1300307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1301307b7a18SHong Zhang mat_a += bs2; 1302307b7a18SHong Zhang } 1303307b7a18SHong Zhang ilen_row++; 1304307b7a18SHong Zhang } 1305307b7a18SHong Zhang for (l=imark; l<nzB; l++) { 1306307b7a18SHong Zhang *mat_j++ = bmap[cworkB[l]]; 1307307b7a18SHong Zhang if (!ijonly) { 1308307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1309307b7a18SHong Zhang mat_a += bs2; 1310307b7a18SHong Zhang } 1311307b7a18SHong Zhang ilen_row++; 1312307b7a18SHong Zhang } 1313307b7a18SHong Zhang } 1314d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1315d5b485abSSatish Balay } 1316d5b485abSSatish Balay } 1317d5b485abSSatish Balay } 1318d5b485abSSatish Balay } 1319d5b485abSSatish Balay 1320d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1321d5b485abSSatish Balay { 1322b24ad042SBarry Smith PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1323b24ad042SBarry Smith PetscInt *imat_j,*imat_i; 13240298fd71SBarry Smith MatScalar *imat_a = NULL,*rbuf4_i = NULL; 1325b24ad042SBarry Smith PetscMPIInt ii; 1326d5b485abSSatish Balay 1327d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 132826fbe8dcSKarl Rupp if (ijonly) ii = tmp2; 132926fbe8dcSKarl Rupp else { 1330b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr); 13317a868f3eSHong Zhang } 1332b24ad042SBarry Smith idex = pa[ii]; 1333999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1334d5b485abSSatish Balay jmax = sbuf1_i[0]; 1335d5b485abSSatish Balay ct1 = 2*jmax + 1; 1336d5b485abSSatish Balay ct2 = 0; 1337b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 1338b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 13397a868f3eSHong Zhang if (!ijonly) rbuf4_i = rbuf4[ii]; 1340d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1341d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 13428fa52d88SHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1343d5b485abSSatish Balay rmap_i = rmap[is_no]; 1344df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[is_no]->data; 1345d5b485abSSatish Balay imat_ilen = mat->ilen; 1346d5b485abSSatish Balay imat_j = mat->j; 1347d5b485abSSatish Balay imat_i = mat->i; 13487a868f3eSHong Zhang if (!ijonly) imat_a = mat->a; 1349d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1350d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1351d5b485abSSatish Balay row = sbuf1_i[ct1]; 1352233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 13538fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 1354233c2ff6SSatish Balay row--; 135526fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1356233c2ff6SSatish Balay #else 1357d5b485abSSatish Balay row = rmap_i[row]; 1358233c2ff6SSatish Balay #endif 1359d5b485abSSatish Balay ilen = imat_ilen[row]; 1360df36731bSBarry Smith mat_i = imat_i[row]; 13617a868f3eSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1362d5b485abSSatish Balay mat_j = imat_j + mat_i; 1363d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1364307b7a18SHong Zhang 1365307b7a18SHong Zhang if (!allcolumns[is_no]) { 1366d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1367233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 13688fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1369233c2ff6SSatish Balay if (tcol) { 1370233c2ff6SSatish Balay #else 1371d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1372233c2ff6SSatish Balay #endif 1373df36731bSBarry Smith *mat_j++ = tcol - 1; 13747a868f3eSHong Zhang if (!ijonly) { 13753eda8832SBarry Smith ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1376549d3d68SSatish Balay mat_a += bs2; 13777a868f3eSHong Zhang } 1378d5b485abSSatish Balay ilen++; 1379d5b485abSSatish Balay } 1380d5b485abSSatish Balay } 1381307b7a18SHong Zhang } else { /* allcolumns */ 1382307b7a18SHong Zhang for (l=0; l<max2; l++,ct2++) { 1383307b7a18SHong Zhang *mat_j++ = rbuf3_i[ct2]; 1384307b7a18SHong Zhang if (!ijonly) { 1385307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1386307b7a18SHong Zhang mat_a += bs2; 1387307b7a18SHong Zhang } 1388307b7a18SHong Zhang ilen++; 1389307b7a18SHong Zhang } 1390307b7a18SHong Zhang } 1391d5b485abSSatish Balay imat_ilen[row] = ilen; 1392d5b485abSSatish Balay } 1393d5b485abSSatish Balay } 1394d5b485abSSatish Balay } 1395d5b485abSSatish Balay } 13967a868f3eSHong Zhang if (!ijonly) { 1397606d414cSSatish Balay ierr = PetscFree(r_status4);CHKERRQ(ierr); 1398606d414cSSatish Balay ierr = PetscFree(r_waits4);CHKERRQ(ierr); 13990c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 1400606d414cSSatish Balay ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1401606d414cSSatish Balay ierr = PetscFree(s_status4);CHKERRQ(ierr); 14027a868f3eSHong Zhang } 1403d5b485abSSatish Balay 1404d5b485abSSatish Balay /* Restore the indices */ 1405d5b485abSSatish Balay for (i=0; i<ismax; i++) { 14064da72fa9SHong Zhang if (!allrows[i]) { 1407d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 14084da72fa9SHong Zhang } 1409307b7a18SHong Zhang if (!allcolumns[i]) { 1410d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1411d5b485abSSatish Balay } 1412307b7a18SHong Zhang } 1413d5b485abSSatish Balay 1414d5b485abSSatish Balay /* Destroy allocated memory */ 141523bdbc58SBarry Smith #if defined(PETSC_USE_CTABLE) 141623bdbc58SBarry Smith ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr); 141723bdbc58SBarry Smith #else 141823bdbc58SBarry Smith ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr); 141923bdbc58SBarry Smith #endif 142023bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 1421606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 1422d5b485abSSatish Balay 142323bdbc58SBarry Smith ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 1424606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 1425606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 1426d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1427606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1428d5b485abSSatish Balay } 1429d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1430606d414cSSatish Balay ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1431d5b485abSSatish Balay } 143223bdbc58SBarry Smith ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr); 1433606d414cSSatish Balay ierr = PetscFree(rbuf3);CHKERRQ(ierr); 1434606d414cSSatish Balay ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1435606d414cSSatish Balay ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 14367a868f3eSHong Zhang if (!ijonly) { 14377a868f3eSHong Zhang for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);} 14387a868f3eSHong Zhang ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1439606d414cSSatish Balay ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1440606d414cSSatish Balay ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 14417a868f3eSHong Zhang } 1442d5b485abSSatish Balay 1443233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1444ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 14458fa52d88SHong Zhang ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr); 1446233c2ff6SSatish Balay } 1447233c2ff6SSatish Balay #endif 14488fa52d88SHong Zhang ierr = PetscFree(rmap);CHKERRQ(ierr); 14498fa52d88SHong Zhang 14508fa52d88SHong Zhang for (i=0; i<ismax; i++) { 14518fa52d88SHong Zhang if (!allcolumns[i]) { 14528fa52d88SHong Zhang #if defined(PETSC_USE_CTABLE) 14538fa52d88SHong Zhang ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr); 14548fa52d88SHong Zhang #else 14558fa52d88SHong Zhang ierr = PetscFree(cmap[i]);CHKERRQ(ierr); 14568fa52d88SHong Zhang #endif 14578fa52d88SHong Zhang } 14588fa52d88SHong Zhang } 14598fa52d88SHong Zhang ierr = PetscFree(cmap);CHKERRQ(ierr); 1460606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1461d5b485abSSatish Balay 1462d5b485abSSatish Balay for (i=0; i<ismax; i++) { 146336f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146436f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1465d5b485abSSatish Balay } 14667a868f3eSHong Zhang 14677a868f3eSHong Zhang c->ijonly = PETSC_FALSE; /* set back to the default */ 14683a40ed3dSBarry Smith PetscFunctionReturn(0); 1469d5b485abSSatish Balay } 1470ca161407SBarry Smith 1471