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 14b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov) 15d5b485abSSatish Balay { 166849ba73SBarry Smith PetscErrorCode ierr; 17d0f46423SBarry Smith PetscInt i,N=C->cmap->N, bs=C->rmap->bs; 1836f4e84dSSatish Balay IS *is_new; 1936f4e84dSSatish Balay 203a40ed3dSBarry Smith PetscFunctionBegin; 21785e854fSJed Brown ierr = PetscMalloc1(imax,&is_new);CHKERRQ(ierr); 22df36731bSBarry Smith /* Convert the indices into block format */ 2305d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr); 2426fbe8dcSKarl Rupp if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n"); 25d5b485abSSatish Balay for (i=0; i<ov; ++i) { 2636f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 27d5b485abSSatish Balay } 286bf464f9SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 2905d8c843SHong Zhang ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr); 306bf464f9SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 31606d414cSSatish Balay ierr = PetscFree(is_new);CHKERRQ(ierr); 323a40ed3dSBarry Smith PetscFunctionReturn(0); 33d5b485abSSatish Balay } 34d5b485abSSatish Balay 35d5b485abSSatish Balay /* 36d5b485abSSatish Balay Sample message format: 37d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 380e9b0e7eSHong Zhang to index sets is[1], is[5] 39d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 40d5b485abSSatish Balay ----------- 41d5b485abSSatish Balay mesg [1] = 1 => is[1] 42d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 43d5b485abSSatish Balay ----------- 44d5b485abSSatish Balay mesg [5] = 5 => is[5] 45d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 46d5b485abSSatish Balay ----------- 47d5b485abSSatish Balay mesg [7] 480e9b0e7eSHong Zhang mesg [n] data(is[1]) 49d5b485abSSatish Balay ----------- 50d5b485abSSatish Balay mesg[n+1] 51d5b485abSSatish Balay mesg[m] data(is[5]) 52d5b485abSSatish Balay ----------- 53d5b485abSSatish Balay 54d5b485abSSatish Balay Notes: 55d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 56d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 57d5b485abSSatish Balay */ 58db41cccfSHong Zhang PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[]) 59d5b485abSSatish Balay { 60df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 615d0c19d7SBarry Smith const PetscInt **idx,*idx_i; 6224c049a4SHong Zhang PetscInt *n,*w3,*w4,**data,len; 636849ba73SBarry Smith PetscErrorCode ierr; 64b24ad042SBarry Smith PetscMPIInt size,rank,tag1,tag2,*w2,*w1,nrqr; 65245d216aSHong Zhang PetscInt Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr; 668f9f447aSHong Zhang PetscInt *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p; 67b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1,*onodes2,*olengths2; 68f1af5d2fSBarry Smith PetscBT *table; 69d5b485abSSatish Balay MPI_Comm comm; 70d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 71d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 728f9f447aSHong Zhang char *t_p; 73d5b485abSSatish Balay 743a40ed3dSBarry Smith PetscFunctionBegin; 75ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 76d5b485abSSatish Balay size = c->size; 77d5b485abSSatish Balay rank = c->rank; 78df36731bSBarry Smith Mbs = c->Mbs; 79d5b485abSSatish Balay 80c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 81c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 82c7dd2924SSatish Balay 83dcca6d9dSJed Brown ierr = PetscMalloc2(imax+1,&idx,imax,&n);CHKERRQ(ierr); 8424c049a4SHong Zhang 85d5b485abSSatish Balay for (i=0; i<imax; i++) { 86d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 87b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 88d5b485abSSatish Balay } 89d5b485abSSatish Balay 90d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 91d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 92884e879aSBarry Smith ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 93d5b485abSSatish Balay for (i=0; i<imax; i++) { 94b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 95d5b485abSSatish Balay idx_i = idx[i]; 96d5b485abSSatish Balay len = n[i]; 97d5b485abSSatish Balay for (j=0; j<len; j++) { 98d5b485abSSatish Balay row = idx_i[j]; 99f23aa3ddSBarry Smith if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 100ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 101d5b485abSSatish Balay w4[proc]++; 102d5b485abSSatish Balay } 103d5b485abSSatish Balay for (j=0; j<size; j++) { 104d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 105d5b485abSSatish Balay } 106d5b485abSSatish Balay } 107d5b485abSSatish Balay 108d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 109d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1100e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 111d5b485abSSatish Balay w3[rank] = 0; 112d5b485abSSatish Balay for (i=0; i<size; i++) { 113d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 114d5b485abSSatish Balay } 115d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 116854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); 117d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 118d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 119d5b485abSSatish Balay } 120d5b485abSSatish Balay 121d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 122d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 123d5b485abSSatish Balay j = pa[i]; 124d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 125d5b485abSSatish Balay msz += w1[j]; 126d5b485abSSatish Balay } 127d5b485abSSatish Balay 128c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 129a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 130c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 131d5b485abSSatish Balay 132c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 133c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 134d5b485abSSatish Balay 135d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 136dcca6d9dSJed Brown ierr = PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr);CHKERRQ(ierr); 13723bdbc58SBarry Smith ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr); 13823bdbc58SBarry Smith ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 139d5b485abSSatish Balay { 140b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 141d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 142d5b485abSSatish Balay j = pa[i]; 143d5b485abSSatish Balay iptr += ict; 144d5b485abSSatish Balay outdat[j] = iptr; 145d5b485abSSatish Balay ict = w1[j]; 146d5b485abSSatish Balay } 147d5b485abSSatish Balay } 148d5b485abSSatish Balay 149d5b485abSSatish Balay /* Form the outgoing messages */ 150d5b485abSSatish Balay /*plug in the headers*/ 151d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 152d5b485abSSatish Balay j = pa[i]; 153d5b485abSSatish Balay outdat[j][0] = 0; 154b24ad042SBarry Smith ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 155d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 156d5b485abSSatish Balay } 157d5b485abSSatish Balay 158d5b485abSSatish Balay /* Memory for doing local proc's work*/ 159d5b485abSSatish Balay { 1601795a4d1SJed Brown ierr = PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mbs*imax,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,&t_p);CHKERRQ(ierr); 161d5b485abSSatish Balay 162bbd702dbSSatish Balay for (i=0; i<imax; i++) { 163b6410449SSatish Balay table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 164bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 165d5b485abSSatish Balay } 166d5b485abSSatish Balay } 167d5b485abSSatish Balay 168d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 169d5b485abSSatish Balay { 170b24ad042SBarry Smith PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 171f1af5d2fSBarry Smith PetscBT table_i; 172d5b485abSSatish Balay 173d5b485abSSatish Balay for (i=0; i<imax; i++) { 174b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 175d5b485abSSatish Balay n_i = n[i]; 176d5b485abSSatish Balay table_i = table[i]; 177d5b485abSSatish Balay idx_i = idx[i]; 178d5b485abSSatish Balay data_i = data[i]; 179d5b485abSSatish Balay isz_i = isz[i]; 180d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 181d5b485abSSatish Balay row = idx_i[j]; 182ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 183d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 184d5b485abSSatish Balay ctr[proc]++; 185d5b485abSSatish Balay *ptr[proc] = row; 186d5b485abSSatish Balay ptr[proc]++; 187d6b45a43SBarry Smith } else { /* Update the local table */ 18826fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 189d5b485abSSatish Balay } 190d5b485abSSatish Balay } 191d5b485abSSatish Balay /* Update the headers for the current IS */ 192d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 193d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 194d5b485abSSatish Balay outdat_j = outdat[j]; 195d5b485abSSatish Balay k = ++outdat_j[0]; 196d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 197d5b485abSSatish Balay outdat_j[2*k-1] = i; 198d5b485abSSatish Balay } 199d5b485abSSatish Balay } 200d5b485abSSatish Balay isz[i] = isz_i; 201d5b485abSSatish Balay } 202d5b485abSSatish Balay } 203d5b485abSSatish Balay 204d5b485abSSatish Balay /* Now post the sends */ 205854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 206d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 207d5b485abSSatish Balay j = pa[i]; 208b24ad042SBarry Smith ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 209d5b485abSSatish Balay } 210d5b485abSSatish Balay 211d5b485abSSatish Balay /* No longer need the original indices*/ 212d5b485abSSatish Balay for (i=0; i<imax; ++i) { 213d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 214d5b485abSSatish Balay } 21524c049a4SHong Zhang ierr = PetscFree2(idx,n);CHKERRQ(ierr); 216d5b485abSSatish Balay 217d5b485abSSatish Balay for (i=0; i<imax; ++i) { 2186bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 219d5b485abSSatish Balay } 220d5b485abSSatish Balay 221d5b485abSSatish Balay /* Do Local work*/ 222df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 223d5b485abSSatish Balay 224d5b485abSSatish Balay /* Receive messages*/ 225854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&recv_status);CHKERRQ(ierr); 2260c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);} 227d5b485abSSatish Balay 228854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_status);CHKERRQ(ierr); 2290c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 230d5b485abSSatish Balay 231d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 23223bdbc58SBarry Smith ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr); 23323bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 234d5b485abSSatish Balay 235854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&xdata);CHKERRQ(ierr); 236854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&isz1);CHKERRQ(ierr); 237df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 238c05d87d6SBarry Smith ierr = PetscFree(rbuf[0]);CHKERRQ(ierr); 239606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 240d5b485abSSatish Balay 241d5b485abSSatish Balay /* Send the data back*/ 242d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 243d5b485abSSatish Balay { 244b24ad042SBarry Smith PetscMPIInt *rw1; 245d5b485abSSatish Balay 246884e879aSBarry Smith ierr = PetscCalloc1(size,&rw1);CHKERRQ(ierr); 247c7dd2924SSatish Balay 248d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 249d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 250e32f2f54SBarry Smith if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 251d5b485abSSatish Balay rw1[proc] = isz1[i]; 252d5b485abSSatish Balay } 253d5b485abSSatish Balay 254c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 255c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 256d5b485abSSatish Balay 257c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 258c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 25903512dc5SSatish Balay ierr = PetscFree(rw1);CHKERRQ(ierr); 260c7dd2924SSatish Balay } 261c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 262c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 263d5b485abSSatish Balay 264d5b485abSSatish Balay /* Now post the sends */ 265854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 266d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 267d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 268b24ad042SBarry Smith ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 269d5b485abSSatish Balay } 270d5b485abSSatish Balay 271d5b485abSSatish Balay /* receive work done on other processors*/ 272d5b485abSSatish Balay { 273b24ad042SBarry Smith PetscMPIInt idex; 274b24ad042SBarry Smith PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 275f1af5d2fSBarry Smith PetscBT table_i; 276d5b485abSSatish Balay MPI_Status *status2; 277d5b485abSSatish Balay 278854ce69bSBarry Smith ierr = PetscMalloc1(PetscMax(nrqr,nrqs)+1,&status2);CHKERRQ(ierr); 279d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 280999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 281d5b485abSSatish Balay /* Process the message*/ 282999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 283d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 284999d9058SBarry Smith jmax = rbuf2[idex][0]; 285d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 286d5b485abSSatish Balay max = rbuf2_i[2*j]; 287d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 288d5b485abSSatish Balay isz_i = isz[is_no]; 289d5b485abSSatish Balay data_i = data[is_no]; 290d5b485abSSatish Balay table_i = table[is_no]; 291d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 292d5b485abSSatish Balay row = rbuf2_i[ct1]; 29326fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 294d5b485abSSatish Balay } 295d5b485abSSatish Balay isz[is_no] = isz_i; 296d5b485abSSatish Balay } 297d5b485abSSatish Balay } 2980c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);} 299606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 300d5b485abSSatish Balay } 301d5b485abSSatish Balay 302d5b485abSSatish Balay for (i=0; i<imax; ++i) { 30370b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 304d5b485abSSatish Balay } 305d5b485abSSatish Balay 306c7dd2924SSatish Balay 307c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 308c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 309c7dd2924SSatish Balay 310606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 311c05d87d6SBarry Smith ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr); 312606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 313606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 314606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 315606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 316606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 3178f9f447aSHong Zhang ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr); 318606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 319606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 320606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 321606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 322606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 3233a40ed3dSBarry Smith PetscFunctionReturn(0); 324d5b485abSSatish Balay } 325d5b485abSSatish Balay 326d5b485abSSatish Balay /* 327df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 328d5b485abSSatish Balay the work on the local processor. 329d5b485abSSatish Balay 330d5b485abSSatish Balay Inputs: 331df36731bSBarry Smith C - MAT_MPIBAIJ; 332d5b485abSSatish Balay imax - total no of index sets processed at a time; 333df36731bSBarry Smith table - an array of char - size = Mbs bits. 334d5b485abSSatish Balay 335d5b485abSSatish Balay Output: 3360e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 337d5b485abSSatish Balay to each index set; 338d5b485abSSatish Balay data - pointer to the solutions 339d5b485abSSatish Balay */ 340b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 341d5b485abSSatish Balay { 342df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 343d5b485abSSatish Balay Mat A = c->A,B = c->B; 344df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 345b24ad042SBarry Smith PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 346b24ad042SBarry Smith PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 347f1af5d2fSBarry Smith PetscBT table_i; 348d5b485abSSatish Balay 3493a40ed3dSBarry Smith PetscFunctionBegin; 350899cda47SBarry Smith rstart = c->rstartbs; 351899cda47SBarry Smith cstart = c->cstartbs; 352d5b485abSSatish Balay ai = a->i; 353df36731bSBarry Smith aj = a->j; 354d5b485abSSatish Balay bi = b->i; 355df36731bSBarry Smith bj = b->j; 356d5b485abSSatish Balay garray = c->garray; 357d5b485abSSatish Balay 358d5b485abSSatish Balay 359d5b485abSSatish Balay for (i=0; i<imax; i++) { 360d5b485abSSatish Balay data_i = data[i]; 361d5b485abSSatish Balay table_i = table[i]; 362d5b485abSSatish Balay isz_i = isz[i]; 363d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 364d5b485abSSatish Balay row = data_i[j] - rstart; 365d5b485abSSatish Balay start = ai[row]; 366d5b485abSSatish Balay end = ai[row+1]; 367d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 368df36731bSBarry Smith val = aj[k] + cstart; 36926fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 370d5b485abSSatish Balay } 371d5b485abSSatish Balay start = bi[row]; 372d5b485abSSatish Balay end = bi[row+1]; 373d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 374df36731bSBarry Smith val = garray[bj[k]]; 37526fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 376d5b485abSSatish Balay } 377d5b485abSSatish Balay } 378d5b485abSSatish Balay isz[i] = isz_i; 379d5b485abSSatish Balay } 3803a40ed3dSBarry Smith PetscFunctionReturn(0); 381d5b485abSSatish Balay } 382d5b485abSSatish Balay /* 383df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 384d5b485abSSatish Balay and return the output 385d5b485abSSatish Balay 386d5b485abSSatish Balay Input: 387d5b485abSSatish Balay C - the matrix 388d5b485abSSatish Balay nrqr - no of messages being processed. 389d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 390d5b485abSSatish Balay 391d5b485abSSatish Balay Output: 392d5b485abSSatish Balay xdata - array of messages to be sent back 393d5b485abSSatish Balay isz1 - size of each message 394d5b485abSSatish Balay 395a8c7a070SBarry Smith For better efficiency perhaps we should malloc separately each xdata[i], 396d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 3970e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of 398d5b485abSSatish Balay memory is used. 399d5b485abSSatish Balay 400d5b485abSSatish Balay */ 401b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 402d5b485abSSatish Balay { 403df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 404d5b485abSSatish Balay Mat A = c->A,B = c->B; 405df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 4066849ba73SBarry Smith PetscErrorCode ierr; 407b24ad042SBarry Smith PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 408b24ad042SBarry Smith PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 409d2a212eaSBarry Smith PetscInt val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 410b24ad042SBarry Smith PetscInt *rbuf_i,kmax,rbuf_0; 411f1af5d2fSBarry Smith PetscBT xtable; 412d5b485abSSatish Balay 4133a40ed3dSBarry Smith PetscFunctionBegin; 414df36731bSBarry Smith Mbs = c->Mbs; 415899cda47SBarry Smith rstart = c->rstartbs; 416899cda47SBarry Smith cstart = c->cstartbs; 417d5b485abSSatish Balay ai = a->i; 418df36731bSBarry Smith aj = a->j; 419d5b485abSSatish Balay bi = b->i; 420df36731bSBarry Smith bj = b->j; 421d5b485abSSatish Balay garray = c->garray; 422d5b485abSSatish Balay 423d5b485abSSatish Balay 424d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 425d5b485abSSatish Balay rbuf_i = rbuf[i]; 426d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 427d5b485abSSatish Balay ct += rbuf_0; 42826fbe8dcSKarl Rupp for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j]; 429d5b485abSSatish Balay } 430d5b485abSSatish Balay 431701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 432701b8814SBarry Smith else max1 = 1; 433d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 434785e854fSJed Brown ierr = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr); 435d5b485abSSatish Balay ++no_malloc; 43653b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr); 437b24ad042SBarry Smith ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr); 438d5b485abSSatish Balay 439d5b485abSSatish Balay ct3 = 0; 440d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 441d5b485abSSatish Balay rbuf_i = rbuf[i]; 442d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 443d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 444d5b485abSSatish Balay ct2 = ct1; 445d5b485abSSatish Balay ct3 += ct1; 446d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 4476831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 448d5b485abSSatish Balay oct2 = ct2; 449d5b485abSSatish Balay kmax = rbuf_i[2*j]; 450d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 451d5b485abSSatish Balay row = rbuf_i[ct1]; 452f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 453d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 454b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 455785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 456b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 457606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 458d5b485abSSatish Balay xdata[0] = tmp; 459d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 46026fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 461d5b485abSSatish Balay } 462d5b485abSSatish Balay xdata[i][ct2++] = row; 463d5b485abSSatish Balay ct3++; 464d5b485abSSatish Balay } 465d5b485abSSatish Balay } 466d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 467d5b485abSSatish Balay row = xdata[i][k] - rstart; 468d5b485abSSatish Balay start = ai[row]; 469d5b485abSSatish Balay end = ai[row+1]; 470d5b485abSSatish Balay for (l=start; l<end; l++) { 471df36731bSBarry Smith val = aj[l] + cstart; 472f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 473d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 474b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 475785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 476b24ad042SBarry Smith ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 477606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 478d5b485abSSatish Balay xdata[0] = tmp; 479d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 48026fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 481d5b485abSSatish Balay } 482d5b485abSSatish Balay xdata[i][ct2++] = val; 483d5b485abSSatish Balay ct3++; 484d5b485abSSatish Balay } 485d5b485abSSatish Balay } 486d5b485abSSatish Balay start = bi[row]; 487d5b485abSSatish Balay end = bi[row+1]; 488d5b485abSSatish Balay for (l=start; l<end; l++) { 489df36731bSBarry Smith val = garray[bj[l]]; 490f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 491d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 492b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 493785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&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 } 505d5b485abSSatish Balay /* Update the header*/ 506d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 507d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 508d5b485abSSatish Balay } 509d5b485abSSatish Balay xdata[i][0] = rbuf_0; 510d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 511d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 512d5b485abSSatish Balay } 51394bacf5dSBarry Smith ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr); 5141e2582c4SBarry Smith ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr); 5153a40ed3dSBarry Smith PetscFunctionReturn(0); 516d5b485abSSatish Balay } 517d5b485abSSatish Balay 518b24ad042SBarry Smith PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 519d5b485abSSatish Balay { 52017df9f7cSHong Zhang IS *isrow_block,*iscol_block; 521cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 5226849ba73SBarry Smith PetscErrorCode ierr; 5234da72fa9SHong Zhang PetscInt nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs; 5244da72fa9SHong Zhang PetscBool colflag,*allcolumns,*allrows; 525a2feddc5SSatish Balay 5263a40ed3dSBarry Smith PetscFunctionBegin; 52717df9f7cSHong Zhang //printf("MatGetSubMatrices_MPIBAIJ scall %d, bs %d\n",scall,bs); 52829dcf524SDmitry Karpeev /* The compression and expansion should be avoided. Doesn't point 52929dcf524SDmitry Karpeev out errors, might change the indices, hence buggey */ 53017df9f7cSHong Zhang ierr = PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block);CHKERRQ(ierr); 53117df9f7cSHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block);CHKERRQ(ierr); 53217df9f7cSHong Zhang ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block);CHKERRQ(ierr); 533cf4f527aSSatish Balay 5345dd0c08cSHong Zhang /* Check for special case: each processor gets entire matrix columns -- rm! */ 535dcca6d9dSJed Brown ierr = PetscMalloc2(ismax+1,&allcolumns,ismax+1,&allrows);CHKERRQ(ierr); 536307b7a18SHong Zhang for (i=0; i<ismax; i++) { 537307b7a18SHong Zhang ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 538307b7a18SHong Zhang ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr); 53926fbe8dcSKarl Rupp if (colflag && ncol == C->cmap->N) allcolumns[i] = PETSC_TRUE; 54026fbe8dcSKarl Rupp else allcolumns[i] = PETSC_FALSE; 5414da72fa9SHong Zhang 5424da72fa9SHong Zhang ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 5434da72fa9SHong Zhang ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr); 54426fbe8dcSKarl Rupp if (colflag && nrow == C->rmap->N) allrows[i] = PETSC_TRUE; 54526fbe8dcSKarl Rupp else allrows[i] = PETSC_FALSE; 546307b7a18SHong Zhang } 547307b7a18SHong Zhang 548cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 54917df9f7cSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 550854ce69bSBarry Smith ierr = PetscMalloc1(ismax+1,submat);CHKERRQ(ierr); 551cf4f527aSSatish Balay } 552cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 553b24ad042SBarry Smith nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 554329f5518SBarry Smith if (!nmax) nmax = 1; 555cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); 556cf4f527aSSatish Balay 557653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 558b2566f29SBarry Smith ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 559cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 560cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 561cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 562cf4f527aSSatish Balay else max_no = ismax-pos; 563a42ab0d6SHong Zhang 564a42ab0d6SHong Zhang PetscBool isbaij; 565a42ab0d6SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)C,MATMPIBAIJ,&isbaij);CHKERRQ(ierr); 566a42ab0d6SHong Zhang if (isbaij) { 56717df9f7cSHong Zhang ierr = MatGetSubMatrices_MPIBAIJ_local_new(C,max_no,isrow_block+pos,iscol_block+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 568a42ab0d6SHong Zhang } else { 569a42ab0d6SHong Zhang ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 570a42ab0d6SHong Zhang } 571cf4f527aSSatish Balay pos += max_no; 572cf4f527aSSatish Balay } 57336f4e84dSSatish Balay 57436f4e84dSSatish Balay for (i=0; i<ismax; i++) { 57517df9f7cSHong Zhang ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr); 57617df9f7cSHong Zhang ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr); 57736f4e84dSSatish Balay } 57817df9f7cSHong Zhang ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr); 5794da72fa9SHong Zhang ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr); 5803a40ed3dSBarry Smith PetscFunctionReturn(0); 581a2feddc5SSatish Balay } 582a2feddc5SSatish Balay 583233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 584e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 585233c2ff6SSatish Balay { 586e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 5874dc2109aSBarry Smith PetscMPIInt fproc; 5884dc2109aSBarry Smith PetscErrorCode ierr; 589233c2ff6SSatish Balay 590233c2ff6SSatish Balay PetscFunctionBegin; 5914dc2109aSBarry Smith ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr); 59223ce1328SBarry Smith if (fproc > size) fproc = size; 593e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 594e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 595233c2ff6SSatish Balay else fproc++; 596233c2ff6SSatish Balay } 597e005ede5SBarry Smith *rank = fproc; 598233c2ff6SSatish Balay PetscFunctionReturn(0); 599233c2ff6SSatish Balay } 600233c2ff6SSatish Balay #endif 601233c2ff6SSatish Balay 602a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 603b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 60417df9f7cSHong Zhang 60517df9f7cSHong Zhang PetscErrorCode MatDestroy_MPIBAIJ_MatGetSubmatrices(Mat C) 60617df9f7cSHong Zhang { 60717df9f7cSHong Zhang PetscErrorCode ierr; 60817df9f7cSHong Zhang Mat_SeqBAIJ *c = (Mat_SeqBAIJ*)C->data; 60917df9f7cSHong Zhang Mat_SubMat *submatj = c->submatis1; 61017df9f7cSHong Zhang PetscInt i; 61117df9f7cSHong Zhang 61217df9f7cSHong Zhang PetscFunctionBegin; 61317df9f7cSHong Zhang if (!submatj->id) { /* delete data that are linked only to submats[id=0] */ 61417df9f7cSHong Zhang ierr = PetscFree4(submatj->sbuf1,submatj->ptr,submatj->tmp,submatj->ctr);CHKERRQ(ierr); 61517df9f7cSHong Zhang 61617df9f7cSHong Zhang for (i=0; i<submatj->nrqr; ++i) { 61717df9f7cSHong Zhang ierr = PetscFree(submatj->sbuf2[i]);CHKERRQ(ierr); 61817df9f7cSHong Zhang } 61917df9f7cSHong Zhang ierr = PetscFree3(submatj->sbuf2,submatj->req_size,submatj->req_source1);CHKERRQ(ierr); 62017df9f7cSHong Zhang 62117df9f7cSHong Zhang if (submatj->rbuf1) { 62217df9f7cSHong Zhang ierr = PetscFree(submatj->rbuf1[0]);CHKERRQ(ierr); 62317df9f7cSHong Zhang ierr = PetscFree(submatj->rbuf1);CHKERRQ(ierr); 62417df9f7cSHong Zhang } 62517df9f7cSHong Zhang 62617df9f7cSHong Zhang for (i=0; i<submatj->nrqs; ++i) { 62717df9f7cSHong Zhang ierr = PetscFree(submatj->rbuf3[i]);CHKERRQ(ierr); 62817df9f7cSHong Zhang } 62917df9f7cSHong Zhang ierr = PetscFree3(submatj->req_source2,submatj->rbuf2,submatj->rbuf3);CHKERRQ(ierr); 63017df9f7cSHong Zhang ierr = PetscFree(submatj->pa);CHKERRQ(ierr); 63117df9f7cSHong Zhang } 63217df9f7cSHong Zhang 63317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 63417df9f7cSHong Zhang ierr = PetscTableDestroy((PetscTable*)&submatj->rmap);CHKERRQ(ierr); 63517df9f7cSHong Zhang if (submatj->cmap_loc) {ierr = PetscFree(submatj->cmap_loc);CHKERRQ(ierr);} 63617df9f7cSHong Zhang ierr = PetscFree(submatj->rmap_loc);CHKERRQ(ierr); 63717df9f7cSHong Zhang #else 63817df9f7cSHong Zhang ierr = PetscFree(submatj->rmap);CHKERRQ(ierr); 63917df9f7cSHong Zhang #endif 64017df9f7cSHong Zhang 64117df9f7cSHong Zhang if (!submatj->allcolumns) { 64217df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 64317df9f7cSHong Zhang ierr = PetscTableDestroy((PetscTable*)&submatj->cmap);CHKERRQ(ierr); 64417df9f7cSHong Zhang #else 64517df9f7cSHong Zhang ierr = PetscFree(submatj->cmap);CHKERRQ(ierr); 64617df9f7cSHong Zhang #endif 64717df9f7cSHong Zhang } 64817df9f7cSHong Zhang ierr = submatj->destroy(C);CHKERRQ(ierr); 64917df9f7cSHong Zhang ierr = PetscFree(submatj->row2proc);CHKERRQ(ierr); 65017df9f7cSHong Zhang 65117df9f7cSHong Zhang ierr = PetscFree(submatj);CHKERRQ(ierr); 65217df9f7cSHong Zhang PetscFunctionReturn(0); 65317df9f7cSHong Zhang } 65417df9f7cSHong Zhang 65517df9f7cSHong Zhang PetscErrorCode MatGetSubMatrices_MPIBAIJ_local_new(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 65617df9f7cSHong Zhang { 65717df9f7cSHong Zhang Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 65817df9f7cSHong Zhang Mat A = c->A; 65917df9f7cSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc; 66017df9f7cSHong Zhang const PetscInt **icol,**irow; 66117df9f7cSHong Zhang PetscInt *nrow,*ncol,start; 66217df9f7cSHong Zhang PetscErrorCode ierr; 66317df9f7cSHong Zhang PetscMPIInt rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr; 66417df9f7cSHong Zhang PetscInt **sbuf1,**sbuf2,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1; 66517df9f7cSHong Zhang PetscInt nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol; 66617df9f7cSHong Zhang PetscInt **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2; 66717df9f7cSHong Zhang PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 66817df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 66917df9f7cSHong Zhang PetscTable *cmap,cmap_i=NULL,*rmap,rmap_i; 67017df9f7cSHong Zhang #else 67117df9f7cSHong Zhang PetscInt **cmap,*cmap_i=NULL,**rmap,*rmap_i; 67217df9f7cSHong Zhang #endif 67317df9f7cSHong Zhang const PetscInt *irow_i; 67417df9f7cSHong Zhang PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 67517df9f7cSHong Zhang MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 67617df9f7cSHong Zhang MPI_Request *r_waits4,*s_waits3,*s_waits4; 67717df9f7cSHong Zhang MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 67817df9f7cSHong Zhang MPI_Status *r_status3,*r_status4,*s_status4; 67917df9f7cSHong Zhang MPI_Comm comm; 68017df9f7cSHong Zhang PetscScalar **rbuf4,*rbuf4_i,**sbuf_aa,*vals,*mat_a,*imat_a,*sbuf_aa_i; 68117df9f7cSHong Zhang PetscMPIInt *onodes1,*olengths1,end; 6829e686edcSHong Zhang PetscInt **row2proc,*row2proc_i,ilen_row,*imat_ilen,*imat_j,*imat_i; 68317df9f7cSHong Zhang Mat_SubMat **smats,*smat_i; 68417df9f7cSHong Zhang PetscBool *issorted,colflag,iscsorted=PETSC_TRUE; 68517df9f7cSHong Zhang PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen; 68617df9f7cSHong Zhang 687a42ab0d6SHong Zhang PetscInt bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs; 688a42ab0d6SHong Zhang PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 689a42ab0d6SHong Zhang PetscInt nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB; 690a42ab0d6SHong Zhang PetscScalar *vworkA,*vworkB,*a_a = a->a,*b_a = b->a; 691a42ab0d6SHong Zhang PetscInt cstart = c->cstartbs,*bmap = c->garray; 692a42ab0d6SHong Zhang 69317df9f7cSHong Zhang PetscFunctionBegin; 69417df9f7cSHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 69517df9f7cSHong Zhang size = c->size; 69617df9f7cSHong Zhang rank = c->rank; 69717df9f7cSHong Zhang 69817df9f7cSHong Zhang ierr = PetscMalloc5(ismax,&smats,ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax,&allcolumns);CHKERRQ(ierr); 69917df9f7cSHong Zhang ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr); 70017df9f7cSHong Zhang 70117df9f7cSHong Zhang for (i=0; i<ismax; i++) { 70217df9f7cSHong Zhang ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr); 70317df9f7cSHong Zhang if (!issorted[i]) iscsorted = issorted[i]; 70417df9f7cSHong Zhang 70517df9f7cSHong Zhang ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr); 70617df9f7cSHong Zhang 70717df9f7cSHong Zhang ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 70817df9f7cSHong Zhang ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 70917df9f7cSHong Zhang 71017df9f7cSHong Zhang /* Check for special case: allcolumn */ 71117df9f7cSHong Zhang ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 71217df9f7cSHong Zhang ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 7135dd0c08cSHong Zhang 714*1abc2f7bSHong Zhang if (colflag && ncol[i] == c->Nbs) { 715*1abc2f7bSHong Zhang allcolumns[i] = PETSC_TRUE; 716*1abc2f7bSHong Zhang icol[i] = NULL; 717*1abc2f7bSHong Zhang printf("[%d] allcolumns[%d] true\n",rank,i); 718*1abc2f7bSHong Zhang } else { 71917df9f7cSHong Zhang allcolumns[i] = PETSC_FALSE; 72017df9f7cSHong Zhang ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 721*1abc2f7bSHong Zhang } 72217df9f7cSHong Zhang } 7235dd0c08cSHong Zhang //printf("[%d] nrow %d, ncol %d\n",rank,nrow[0],ncol[0]); 72417df9f7cSHong Zhang 72517df9f7cSHong Zhang if (scall == MAT_REUSE_MATRIX) { 72617df9f7cSHong Zhang /* Assumes new rows are same length as the old rows */ 72717df9f7cSHong Zhang for (i=0; i<ismax; i++) { 72817df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)(submats[i]->data); 7299e686edcSHong Zhang if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 73017df9f7cSHong Zhang 73117df9f7cSHong Zhang /* Initial matrix as if empty */ 7329e686edcSHong Zhang ierr = PetscMemzero(subc->ilen,subc->mbs*sizeof(PetscInt));CHKERRQ(ierr); 73317df9f7cSHong Zhang 73417df9f7cSHong Zhang /* Initial matrix as if empty */ 73517df9f7cSHong Zhang submats[i]->factortype = C->factortype; 73617df9f7cSHong Zhang 73717df9f7cSHong Zhang smat_i = subc->submatis1; 73817df9f7cSHong Zhang smats[i] = smat_i; 73917df9f7cSHong Zhang 74017df9f7cSHong Zhang nrqs = smat_i->nrqs; 74117df9f7cSHong Zhang nrqr = smat_i->nrqr; 74217df9f7cSHong Zhang rbuf1 = smat_i->rbuf1; 74317df9f7cSHong Zhang rbuf2 = smat_i->rbuf2; 74417df9f7cSHong Zhang rbuf3 = smat_i->rbuf3; 74517df9f7cSHong Zhang req_source2 = smat_i->req_source2; 74617df9f7cSHong Zhang 74717df9f7cSHong Zhang sbuf1 = smat_i->sbuf1; 74817df9f7cSHong Zhang sbuf2 = smat_i->sbuf2; 74917df9f7cSHong Zhang ptr = smat_i->ptr; 75017df9f7cSHong Zhang tmp = smat_i->tmp; 75117df9f7cSHong Zhang ctr = smat_i->ctr; 75217df9f7cSHong Zhang 75317df9f7cSHong Zhang pa = smat_i->pa; 75417df9f7cSHong Zhang req_size = smat_i->req_size; 75517df9f7cSHong Zhang req_source1 = smat_i->req_source1; 75617df9f7cSHong Zhang 75717df9f7cSHong Zhang allcolumns[i] = smat_i->allcolumns; 75817df9f7cSHong Zhang row2proc[i] = smat_i->row2proc; 75917df9f7cSHong Zhang rmap[i] = smat_i->rmap; 76017df9f7cSHong Zhang cmap[i] = smat_i->cmap; 76117df9f7cSHong Zhang } 76217df9f7cSHong Zhang } else { /* scall == MAT_INITIAL_MATRIX */ 76317df9f7cSHong Zhang /* Get some new tags to keep the communication clean */ 76417df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 76517df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 76617df9f7cSHong Zhang 76717df9f7cSHong Zhang /* evaluate communication - mesg to who, length of mesg, and buffer space 76817df9f7cSHong Zhang required. Based on this, buffers are allocated, and data copied into them*/ 76917df9f7cSHong Zhang ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); /* mesg size, initialize work vectors */ 77017df9f7cSHong Zhang 77117df9f7cSHong Zhang for (i=0; i<ismax; i++) { 77217df9f7cSHong Zhang jmax = nrow[i]; 77317df9f7cSHong Zhang irow_i = irow[i]; 77417df9f7cSHong Zhang 77517df9f7cSHong Zhang ierr = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr); 77617df9f7cSHong Zhang row2proc[i] = row2proc_i; 77717df9f7cSHong Zhang 77817df9f7cSHong Zhang if (issorted[i]) proc = 0; 77917df9f7cSHong Zhang for (j=0; j<jmax; j++) { 78017df9f7cSHong Zhang if (!issorted[i]) proc = 0; 78117df9f7cSHong Zhang row = irow_i[j]; 782a42ab0d6SHong Zhang while (row >= c->rangebs[proc+1]) proc++; 78317df9f7cSHong Zhang w4[proc]++; 78417df9f7cSHong Zhang row2proc_i[j] = proc; /* map row index to proc */ 78517df9f7cSHong Zhang } 78617df9f7cSHong Zhang for (j=0; j<size; j++) { 78717df9f7cSHong Zhang if (w4[j]) { w1[j] += w4[j]; w3[j]++; w4[j] = 0;} 78817df9f7cSHong Zhang } 78917df9f7cSHong Zhang } 79017df9f7cSHong Zhang 79117df9f7cSHong Zhang nrqs = 0; /* no of outgoing messages */ 79217df9f7cSHong Zhang msz = 0; /* total mesg length (for all procs) */ 79317df9f7cSHong Zhang w1[rank] = 0; /* no mesg sent to self */ 79417df9f7cSHong Zhang w3[rank] = 0; 79517df9f7cSHong Zhang for (i=0; i<size; i++) { 79617df9f7cSHong Zhang if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 79717df9f7cSHong Zhang } 79817df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 79917df9f7cSHong Zhang for (i=0,j=0; i<size; i++) { 80017df9f7cSHong Zhang if (w1[i]) { pa[j] = i; j++; } 80117df9f7cSHong Zhang } 80217df9f7cSHong Zhang 80317df9f7cSHong Zhang /* Each message would have a header = 1 + 2*(no of IS) + data */ 80417df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 80517df9f7cSHong Zhang j = pa[i]; 80617df9f7cSHong Zhang w1[j] += w2[j] + 2* w3[j]; 80717df9f7cSHong Zhang msz += w1[j]; 80817df9f7cSHong Zhang } 80917df9f7cSHong Zhang ierr = PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);CHKERRQ(ierr); 81017df9f7cSHong Zhang 81117df9f7cSHong Zhang /* Determine the number of messages to expect, their lengths, from from-ids */ 81217df9f7cSHong Zhang ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 81317df9f7cSHong Zhang ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 81417df9f7cSHong Zhang 81517df9f7cSHong Zhang /* Now post the Irecvs corresponding to these messages */ 81617df9f7cSHong Zhang tag0 = ((PetscObject)C)->tag; 81717df9f7cSHong Zhang ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 8185dd0c08cSHong Zhang /* printf("[%d] nrqs %d, nrqr %d\n",rank,nrqs,nrqr); */ 81917df9f7cSHong Zhang 82017df9f7cSHong Zhang ierr = PetscFree(onodes1);CHKERRQ(ierr); 82117df9f7cSHong Zhang ierr = PetscFree(olengths1);CHKERRQ(ierr); 82217df9f7cSHong Zhang 82317df9f7cSHong Zhang /* Allocate Memory for outgoing messages */ 82417df9f7cSHong Zhang ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 82517df9f7cSHong Zhang ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 82617df9f7cSHong Zhang ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 82717df9f7cSHong Zhang 82817df9f7cSHong Zhang { 82917df9f7cSHong Zhang PetscInt *iptr = tmp; 83017df9f7cSHong Zhang k = 0; 83117df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 83217df9f7cSHong Zhang j = pa[i]; 83317df9f7cSHong Zhang iptr += k; 83417df9f7cSHong Zhang sbuf1[j] = iptr; 83517df9f7cSHong Zhang k = w1[j]; 83617df9f7cSHong Zhang } 83717df9f7cSHong Zhang } 83817df9f7cSHong Zhang 83917df9f7cSHong Zhang /* Form the outgoing messages. Initialize the header space */ 84017df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 84117df9f7cSHong Zhang j = pa[i]; 84217df9f7cSHong Zhang sbuf1[j][0] = 0; 84317df9f7cSHong Zhang ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 84417df9f7cSHong Zhang ptr[j] = sbuf1[j] + 2*w3[j] + 1; 84517df9f7cSHong Zhang } 84617df9f7cSHong Zhang 84717df9f7cSHong Zhang /* Parse the isrow and copy data into outbuf */ 84817df9f7cSHong Zhang for (i=0; i<ismax; i++) { 84917df9f7cSHong Zhang row2proc_i = row2proc[i]; 85017df9f7cSHong Zhang ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 85117df9f7cSHong Zhang irow_i = irow[i]; 85217df9f7cSHong Zhang jmax = nrow[i]; 85317df9f7cSHong Zhang for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 85417df9f7cSHong Zhang proc = row2proc_i[j]; 85517df9f7cSHong Zhang if (proc != rank) { /* copy to the outgoing buf*/ 85617df9f7cSHong Zhang ctr[proc]++; 85717df9f7cSHong Zhang *ptr[proc] = irow_i[j]; 85817df9f7cSHong Zhang ptr[proc]++; 85917df9f7cSHong Zhang } 86017df9f7cSHong Zhang } 86117df9f7cSHong Zhang /* Update the headers for the current IS */ 86217df9f7cSHong Zhang for (j=0; j<size; j++) { /* Can Optimise this loop too */ 86317df9f7cSHong Zhang if ((ctr_j = ctr[j])) { 86417df9f7cSHong Zhang sbuf1_j = sbuf1[j]; 86517df9f7cSHong Zhang k = ++sbuf1_j[0]; 86617df9f7cSHong Zhang sbuf1_j[2*k] = ctr_j; 86717df9f7cSHong Zhang sbuf1_j[2*k-1] = i; 86817df9f7cSHong Zhang } 86917df9f7cSHong Zhang } 87017df9f7cSHong Zhang } 87117df9f7cSHong Zhang 87217df9f7cSHong Zhang /* Now post the sends */ 87317df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 87417df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 87517df9f7cSHong Zhang j = pa[i]; 87617df9f7cSHong Zhang ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 87717df9f7cSHong Zhang } 87817df9f7cSHong Zhang 87917df9f7cSHong Zhang /* Post Receives to capture the buffer size */ 88017df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 88117df9f7cSHong Zhang ierr = PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);CHKERRQ(ierr); 88217df9f7cSHong Zhang rbuf2[0] = tmp + msz; 88317df9f7cSHong Zhang for (i=1; i<nrqs; ++i) { 88417df9f7cSHong Zhang rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 88517df9f7cSHong Zhang } 88617df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 88717df9f7cSHong Zhang j = pa[i]; 88817df9f7cSHong Zhang ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRQ(ierr); 88917df9f7cSHong Zhang } 89017df9f7cSHong Zhang 89117df9f7cSHong Zhang /* Send to other procs the buf size they should allocate */ 89217df9f7cSHong Zhang /* Receive messages*/ 89317df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 89417df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 89517df9f7cSHong Zhang ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr); 89617df9f7cSHong Zhang { 897a42ab0d6SHong Zhang PetscInt *sAi = a->i,*sBi = b->i,id; 89817df9f7cSHong Zhang PetscInt *sbuf2_i; 89917df9f7cSHong Zhang 90017df9f7cSHong Zhang ierr = MPI_Waitall(nrqr,r_waits1,r_status1);CHKERRQ(ierr); 90117df9f7cSHong Zhang for (i=0; i<nrqr; ++i) { 90217df9f7cSHong Zhang req_size[i] = 0; 90317df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 90417df9f7cSHong Zhang start = 2*rbuf1_i[0] + 1; 90517df9f7cSHong Zhang ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 90617df9f7cSHong Zhang ierr = PetscMalloc1(end+1,&sbuf2[i]);CHKERRQ(ierr); 90717df9f7cSHong Zhang sbuf2_i = sbuf2[i]; 90817df9f7cSHong Zhang for (j=start; j<end; j++) { 90917df9f7cSHong Zhang id = rbuf1_i[j] - rstart; 91017df9f7cSHong Zhang ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 91117df9f7cSHong Zhang sbuf2_i[j] = ncols; 91217df9f7cSHong Zhang req_size[i] += ncols; 91317df9f7cSHong Zhang } 91417df9f7cSHong Zhang req_source1[i] = r_status1[i].MPI_SOURCE; 91517df9f7cSHong Zhang /* form the header */ 91617df9f7cSHong Zhang sbuf2_i[0] = req_size[i]; 91717df9f7cSHong Zhang for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 91817df9f7cSHong Zhang 91917df9f7cSHong Zhang ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRQ(ierr); 92017df9f7cSHong Zhang } 92117df9f7cSHong Zhang } 92217df9f7cSHong Zhang ierr = PetscFree(r_status1);CHKERRQ(ierr); 92317df9f7cSHong Zhang ierr = PetscFree(r_waits1);CHKERRQ(ierr); 92417df9f7cSHong Zhang ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 92517df9f7cSHong Zhang 92617df9f7cSHong Zhang /* Receive messages*/ 92717df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 92817df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 92917df9f7cSHong Zhang 93017df9f7cSHong Zhang ierr = MPI_Waitall(nrqs,r_waits2,r_status2);CHKERRQ(ierr); 93117df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 93217df9f7cSHong Zhang ierr = PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);CHKERRQ(ierr); 93317df9f7cSHong Zhang req_source2[i] = r_status2[i].MPI_SOURCE; 93417df9f7cSHong Zhang ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRQ(ierr); 93517df9f7cSHong Zhang } 93617df9f7cSHong Zhang ierr = PetscFree(r_status2);CHKERRQ(ierr); 93717df9f7cSHong Zhang ierr = PetscFree(r_waits2);CHKERRQ(ierr); 93817df9f7cSHong Zhang 93917df9f7cSHong Zhang /* Wait on sends1 and sends2 */ 94017df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 94117df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 94217df9f7cSHong Zhang 94317df9f7cSHong Zhang if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 94417df9f7cSHong Zhang if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 94517df9f7cSHong Zhang ierr = PetscFree(s_status1);CHKERRQ(ierr); 94617df9f7cSHong Zhang ierr = PetscFree(s_status2);CHKERRQ(ierr); 94717df9f7cSHong Zhang ierr = PetscFree(s_waits1);CHKERRQ(ierr); 94817df9f7cSHong Zhang ierr = PetscFree(s_waits2);CHKERRQ(ierr); 94917df9f7cSHong Zhang 95017df9f7cSHong Zhang /* Now allocate sending buffers for a->j, and send them off */ 95117df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 95217df9f7cSHong Zhang for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 95317df9f7cSHong Zhang ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 95417df9f7cSHong Zhang for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 95517df9f7cSHong Zhang 95617df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 95717df9f7cSHong Zhang { 95817df9f7cSHong Zhang 95917df9f7cSHong Zhang for (i=0; i<nrqr; i++) { 96017df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 96117df9f7cSHong Zhang sbuf_aj_i = sbuf_aj[i]; 96217df9f7cSHong Zhang ct1 = 2*rbuf1_i[0] + 1; 96317df9f7cSHong Zhang ct2 = 0; 96417df9f7cSHong Zhang for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 96517df9f7cSHong Zhang kmax = rbuf1[i][2*j]; 96617df9f7cSHong Zhang for (k=0; k<kmax; k++,ct1++) { 96717df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 96817df9f7cSHong Zhang nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 96917df9f7cSHong Zhang ncols = nzA + nzB; 97017df9f7cSHong Zhang cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 97117df9f7cSHong Zhang 97217df9f7cSHong Zhang /* load the column indices for this row into cols */ 97317df9f7cSHong Zhang cols = sbuf_aj_i + ct2; 97417df9f7cSHong Zhang for (l=0; l<nzB; l++) { 975a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 976a42ab0d6SHong Zhang else break; 97717df9f7cSHong Zhang } 978a42ab0d6SHong Zhang imark = l; 979a42ab0d6SHong Zhang for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];} 980a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 98117df9f7cSHong Zhang ct2 += ncols; 98217df9f7cSHong Zhang } 98317df9f7cSHong Zhang } 98417df9f7cSHong Zhang ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRQ(ierr); 98517df9f7cSHong Zhang } 98617df9f7cSHong Zhang } 98717df9f7cSHong Zhang ierr = PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);CHKERRQ(ierr); 98817df9f7cSHong Zhang 98917df9f7cSHong Zhang /* create col map: global col of C -> local col of submatrices */ 99017df9f7cSHong Zhang const PetscInt *icol_i; 99117df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 99217df9f7cSHong Zhang for (i=0; i<ismax; i++) { 99317df9f7cSHong Zhang if (!allcolumns[i]) { 994a42ab0d6SHong Zhang ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 99517df9f7cSHong Zhang 99617df9f7cSHong Zhang jmax = ncol[i]; 99717df9f7cSHong Zhang icol_i = icol[i]; 99817df9f7cSHong Zhang cmap_i = cmap[i]; 99917df9f7cSHong Zhang for (j=0; j<jmax; j++) { 100017df9f7cSHong Zhang ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 100117df9f7cSHong Zhang } 100217df9f7cSHong Zhang } else cmap[i] = NULL; 100317df9f7cSHong Zhang } 100417df9f7cSHong Zhang #else 100517df9f7cSHong Zhang for (i=0; i<ismax; i++) { 100617df9f7cSHong Zhang if (!allcolumns[i]) { 1007a42ab0d6SHong Zhang ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 100817df9f7cSHong Zhang jmax = ncol[i]; 100917df9f7cSHong Zhang icol_i = icol[i]; 101017df9f7cSHong Zhang cmap_i = cmap[i]; 1011a42ab0d6SHong Zhang for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1; 101217df9f7cSHong Zhang } else cmap[i] = NULL; 101317df9f7cSHong Zhang } 101417df9f7cSHong Zhang #endif 101517df9f7cSHong Zhang 101617df9f7cSHong Zhang /* Create lens which is required for MatCreate... */ 101717df9f7cSHong Zhang for (i=0,j=0; i<ismax; i++) j += nrow[i]; 101817df9f7cSHong Zhang ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr); 101917df9f7cSHong Zhang 102017df9f7cSHong Zhang if (ismax) { 102117df9f7cSHong Zhang ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr); 102217df9f7cSHong Zhang } 102317df9f7cSHong Zhang for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 102417df9f7cSHong Zhang 102517df9f7cSHong Zhang /* Update lens from local data */ 102617df9f7cSHong Zhang for (i=0; i<ismax; i++) { 102717df9f7cSHong Zhang row2proc_i = row2proc[i]; 102817df9f7cSHong Zhang jmax = nrow[i]; 102917df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 103017df9f7cSHong Zhang irow_i = irow[i]; 103117df9f7cSHong Zhang lens_i = lens[i]; 103217df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1033a42ab0d6SHong Zhang row = irow_i[j]; /* global blocked row of C */ 103417df9f7cSHong Zhang proc = row2proc_i[j]; 103517df9f7cSHong Zhang if (proc == rank) { 1036a42ab0d6SHong Zhang /* Get indices from matA and then from matB */ 103717df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1038a42ab0d6SHong Zhang PetscInt tt; 1039a42ab0d6SHong Zhang #endif 1040a42ab0d6SHong Zhang row = row - rstart; 1041a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1042a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 1043a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 1044a42ab0d6SHong Zhang cworkB = b_j + b_i[row]; 1045a42ab0d6SHong Zhang 1046a42ab0d6SHong Zhang if (!allcolumns[i]) { 1047a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1048a42ab0d6SHong Zhang for (k=0; k<nzA; k++) { 1049a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1050a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1051a42ab0d6SHong Zhang } 1052a42ab0d6SHong Zhang for (k=0; k<nzB; k++) { 1053a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1054a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1055a42ab0d6SHong Zhang } 1056a42ab0d6SHong Zhang 1057a42ab0d6SHong Zhang #else 1058a42ab0d6SHong Zhang for (k=0; k<nzA; k++) { 1059a42ab0d6SHong Zhang if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1060a42ab0d6SHong Zhang } 1061a42ab0d6SHong Zhang for (k=0; k<nzB; k++) { 1062a42ab0d6SHong Zhang if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1063a42ab0d6SHong Zhang } 1064a42ab0d6SHong Zhang #endif 1065a42ab0d6SHong Zhang } else { /* allcolumns */ 1066a42ab0d6SHong Zhang lens_i[j] = nzA + nzB; 1067a42ab0d6SHong Zhang } 106817df9f7cSHong Zhang } 106917df9f7cSHong Zhang } 107017df9f7cSHong Zhang } 107117df9f7cSHong Zhang 107217df9f7cSHong Zhang /* Create row map: global row of C -> local row of submatrices */ 107317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 107417df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1075a42ab0d6SHong Zhang ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 107617df9f7cSHong Zhang irow_i = irow[i]; 107717df9f7cSHong Zhang jmax = nrow[i]; 107817df9f7cSHong Zhang for (j=0; j<jmax; j++) { 107917df9f7cSHong Zhang ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 108017df9f7cSHong Zhang } 108117df9f7cSHong Zhang } 108217df9f7cSHong Zhang #else 108317df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1084a42ab0d6SHong Zhang ierr = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr); 108517df9f7cSHong Zhang rmap_i = rmap[i]; 108617df9f7cSHong Zhang irow_i = irow[i]; 108717df9f7cSHong Zhang jmax = nrow[i]; 108817df9f7cSHong Zhang for (j=0; j<jmax; j++) { 108917df9f7cSHong Zhang rmap_i[irow_i[j]] = j; 109017df9f7cSHong Zhang } 109117df9f7cSHong Zhang } 109217df9f7cSHong Zhang #endif 109317df9f7cSHong Zhang 109417df9f7cSHong Zhang /* Update lens from offproc data */ 109517df9f7cSHong Zhang { 109617df9f7cSHong Zhang PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 109717df9f7cSHong Zhang 109817df9f7cSHong Zhang ierr = MPI_Waitall(nrqs,r_waits3,r_status3);CHKERRQ(ierr); 109917df9f7cSHong Zhang for (tmp2=0; tmp2<nrqs; tmp2++) { 110017df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 110117df9f7cSHong Zhang jmax = sbuf1_i[0]; 110217df9f7cSHong Zhang ct1 = 2*jmax+1; 110317df9f7cSHong Zhang ct2 = 0; 110417df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 110517df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 110617df9f7cSHong Zhang for (j=1; j<=jmax; j++) { 110717df9f7cSHong Zhang is_no = sbuf1_i[2*j-1]; 110817df9f7cSHong Zhang max1 = sbuf1_i[2*j]; 110917df9f7cSHong Zhang lens_i = lens[is_no]; 111017df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 111117df9f7cSHong Zhang rmap_i = rmap[is_no]; 111217df9f7cSHong Zhang for (k=0; k<max1; k++,ct1++) { 111317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 111417df9f7cSHong Zhang ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 111517df9f7cSHong Zhang row--; 111617df9f7cSHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 111717df9f7cSHong Zhang #else 111817df9f7cSHong Zhang row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 111917df9f7cSHong Zhang #endif 112017df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 112117df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 112217df9f7cSHong Zhang if (!allcolumns[is_no]) { 112317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 112417df9f7cSHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 112517df9f7cSHong Zhang #else 112617df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 112717df9f7cSHong Zhang #endif 112817df9f7cSHong Zhang if (tcol) lens_i[row]++; 112917df9f7cSHong Zhang } else { /* allcolumns */ 113017df9f7cSHong Zhang lens_i[row]++; /* lens_i[row] += max2 ? */ 113117df9f7cSHong Zhang } 113217df9f7cSHong Zhang } 113317df9f7cSHong Zhang } 113417df9f7cSHong Zhang } 113517df9f7cSHong Zhang } 113617df9f7cSHong Zhang } 113717df9f7cSHong Zhang ierr = PetscFree(r_waits3);CHKERRQ(ierr); 113817df9f7cSHong Zhang if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 113917df9f7cSHong Zhang ierr = PetscFree2(r_status3,s_status3);CHKERRQ(ierr); 114017df9f7cSHong Zhang ierr = PetscFree(s_waits3);CHKERRQ(ierr); 114117df9f7cSHong Zhang 114217df9f7cSHong Zhang /* Create the submatrices */ 114317df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1144a42ab0d6SHong Zhang PetscInt bs_tmp; 1145a42ab0d6SHong Zhang if (ijonly) bs_tmp = 1; 1146a42ab0d6SHong Zhang else bs_tmp = bs; 114717df9f7cSHong Zhang 114817df9f7cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 1149a42ab0d6SHong Zhang ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 115017df9f7cSHong Zhang 115117df9f7cSHong Zhang ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 1152a42ab0d6SHong Zhang ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 1153a42ab0d6SHong Zhang ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 115417df9f7cSHong Zhang 115517df9f7cSHong Zhang /* create struct Mat_SubMat and attached it to submat */ 115617df9f7cSHong Zhang ierr = PetscNew(&smat_i);CHKERRQ(ierr); 115717df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 115817df9f7cSHong Zhang subc->submatis1 = smat_i; 115917df9f7cSHong Zhang smats[i] = smat_i; 116017df9f7cSHong Zhang 116117df9f7cSHong Zhang smat_i->destroy = submats[i]->ops->destroy; 116217df9f7cSHong Zhang submats[i]->ops->destroy = MatDestroy_MPIBAIJ_MatGetSubmatrices; 116317df9f7cSHong Zhang submats[i]->factortype = C->factortype; 116417df9f7cSHong Zhang 116517df9f7cSHong Zhang smat_i->id = i; 116617df9f7cSHong Zhang smat_i->nrqs = nrqs; 116717df9f7cSHong Zhang smat_i->nrqr = nrqr; 116817df9f7cSHong Zhang smat_i->rbuf1 = rbuf1; 116917df9f7cSHong Zhang smat_i->rbuf2 = rbuf2; 117017df9f7cSHong Zhang smat_i->rbuf3 = rbuf3; 117117df9f7cSHong Zhang smat_i->sbuf2 = sbuf2; 117217df9f7cSHong Zhang smat_i->req_source2 = req_source2; 117317df9f7cSHong Zhang 117417df9f7cSHong Zhang smat_i->sbuf1 = sbuf1; 117517df9f7cSHong Zhang smat_i->ptr = ptr; 117617df9f7cSHong Zhang smat_i->tmp = tmp; 117717df9f7cSHong Zhang smat_i->ctr = ctr; 117817df9f7cSHong Zhang 117917df9f7cSHong Zhang smat_i->pa = pa; 118017df9f7cSHong Zhang smat_i->req_size = req_size; 118117df9f7cSHong Zhang smat_i->req_source1 = req_source1; 118217df9f7cSHong Zhang 118317df9f7cSHong Zhang smat_i->allcolumns = allcolumns[i]; 118417df9f7cSHong Zhang smat_i->singleis = PETSC_FALSE; 118517df9f7cSHong Zhang smat_i->row2proc = row2proc[i]; 118617df9f7cSHong Zhang smat_i->rmap = rmap[i]; 118717df9f7cSHong Zhang smat_i->cmap = cmap[i]; 118817df9f7cSHong Zhang } 118917df9f7cSHong Zhang 119017df9f7cSHong Zhang if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);} 119117df9f7cSHong Zhang ierr = PetscFree(lens);CHKERRQ(ierr); 119217df9f7cSHong Zhang ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 119317df9f7cSHong Zhang ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 119417df9f7cSHong Zhang 119517df9f7cSHong Zhang } /* endof scall == MAT_INITIAL_MATRIX */ 119617df9f7cSHong Zhang 119717df9f7cSHong Zhang /* Post recv matrix values */ 119817df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr); 119917df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 120017df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 120117df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 120217df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 120317df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 1204a42ab0d6SHong Zhang ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr); 1205a42ab0d6SHong Zhang ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRQ(ierr); 120617df9f7cSHong Zhang } 120717df9f7cSHong Zhang 120817df9f7cSHong Zhang /* Allocate sending buffers for a->a, and send them off */ 120917df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 121017df9f7cSHong Zhang for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 12119e686edcSHong Zhang 1212a42ab0d6SHong Zhang ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 1213a42ab0d6SHong Zhang for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 121417df9f7cSHong Zhang 121517df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 121617df9f7cSHong Zhang 121717df9f7cSHong Zhang for (i=0; i<nrqr; i++) { 121817df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 121917df9f7cSHong Zhang sbuf_aa_i = sbuf_aa[i]; 122017df9f7cSHong Zhang ct1 = 2*rbuf1_i[0]+1; 122117df9f7cSHong Zhang ct2 = 0; 122217df9f7cSHong Zhang for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 122317df9f7cSHong Zhang kmax = rbuf1_i[2*j]; 122417df9f7cSHong Zhang for (k=0; k<kmax; k++,ct1++) { 122517df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 1226a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1227a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 122817df9f7cSHong Zhang ncols = nzA + nzB; 122917df9f7cSHong Zhang cworkB = b_j + b_i[row]; 1230a42ab0d6SHong Zhang vworkA = a_a + a_i[row]*bs2; 1231a42ab0d6SHong Zhang vworkB = b_a + b_i[row]*bs2; 123217df9f7cSHong Zhang 123317df9f7cSHong Zhang /* load the column values for this row into vals*/ 1234a42ab0d6SHong Zhang vals = sbuf_aa_i+ct2*bs2; 1235a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1236a42ab0d6SHong Zhang if ((bmap[cworkB[l]]) < cstart) { 1237a42ab0d6SHong Zhang ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1238a42ab0d6SHong Zhang } else break; 1239a42ab0d6SHong Zhang } 1240a42ab0d6SHong Zhang imark = l; 1241a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1242a42ab0d6SHong Zhang ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);//error in proc[1]??? 1243a42ab0d6SHong Zhang } 1244a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1245a42ab0d6SHong Zhang ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1246a42ab0d6SHong Zhang } 124717df9f7cSHong Zhang 124817df9f7cSHong Zhang ct2 += ncols; 124917df9f7cSHong Zhang } 125017df9f7cSHong Zhang } 1251a42ab0d6SHong Zhang ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRQ(ierr); 125217df9f7cSHong Zhang } 125317df9f7cSHong Zhang 125417df9f7cSHong Zhang if (!ismax) { 125517df9f7cSHong Zhang ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 125617df9f7cSHong Zhang ierr = PetscFree(rbuf1);CHKERRQ(ierr); 125717df9f7cSHong Zhang } 125817df9f7cSHong Zhang 125917df9f7cSHong Zhang /* Assemble the matrices */ 126017df9f7cSHong Zhang /* First assemble the local rows */ 126117df9f7cSHong Zhang for (i=0; i<ismax; i++) { 126217df9f7cSHong Zhang row2proc_i = row2proc[i]; 126317df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 126417df9f7cSHong Zhang imat_ilen = subc->ilen; 126517df9f7cSHong Zhang imat_j = subc->j; 126617df9f7cSHong Zhang imat_i = subc->i; 126717df9f7cSHong Zhang imat_a = subc->a; 126817df9f7cSHong Zhang 126917df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 127017df9f7cSHong Zhang rmap_i = rmap[i]; 127117df9f7cSHong Zhang irow_i = irow[i]; 127217df9f7cSHong Zhang jmax = nrow[i]; 127317df9f7cSHong Zhang for (j=0; j<jmax; j++) { 127417df9f7cSHong Zhang row = irow_i[j]; 127517df9f7cSHong Zhang proc = row2proc_i[j]; 1276a42ab0d6SHong Zhang 127717df9f7cSHong Zhang if (proc == rank) { 1278a42ab0d6SHong Zhang 1279a42ab0d6SHong Zhang row = row - rstart; 1280a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1281a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 1282a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 1283a42ab0d6SHong Zhang cworkB = b_j + b_i[row]; 1284a42ab0d6SHong Zhang if (!ijonly) { 1285a42ab0d6SHong Zhang vworkA = a_a + a_i[row]*bs2; 1286a42ab0d6SHong Zhang vworkB = b_a + b_i[row]*bs2; 1287a42ab0d6SHong Zhang } 1288a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1289a42ab0d6SHong Zhang ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1290a42ab0d6SHong Zhang row--; 1291a42ab0d6SHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1292a42ab0d6SHong Zhang #else 1293a42ab0d6SHong Zhang row = rmap_i[row + rstart]; 1294a42ab0d6SHong Zhang #endif 1295a42ab0d6SHong Zhang mat_i = imat_i[row]; 1296a42ab0d6SHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1297a42ab0d6SHong Zhang mat_j = imat_j + mat_i; 1298a42ab0d6SHong Zhang ilen_row = imat_ilen[row]; 1299a42ab0d6SHong Zhang 1300a42ab0d6SHong Zhang /* load the column indices for this row into cols*/ 1301a42ab0d6SHong Zhang if (!allcolumns[i]) { 1302a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1303a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1304a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1305a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1306a42ab0d6SHong Zhang if (tcol) { 1307a42ab0d6SHong Zhang #else 1308a42ab0d6SHong Zhang if ((tcol = cmap_i[ctmp])) { 1309a42ab0d6SHong Zhang #endif 1310a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1311a42ab0d6SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1312a42ab0d6SHong Zhang mat_a += bs2; 1313a42ab0d6SHong Zhang ilen_row++; 1314a42ab0d6SHong Zhang } 1315a42ab0d6SHong Zhang } else break; 1316a42ab0d6SHong Zhang } 1317a42ab0d6SHong Zhang imark = l; 1318a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1319a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1320a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1321a42ab0d6SHong Zhang if (tcol) { 1322a42ab0d6SHong Zhang #else 1323a42ab0d6SHong Zhang if ((tcol = cmap_i[cstart + cworkA[l]])) { 1324a42ab0d6SHong Zhang #endif 1325a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1326a42ab0d6SHong Zhang if (!ijonly) { 1327a42ab0d6SHong Zhang ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1328a42ab0d6SHong Zhang mat_a += bs2; 1329a42ab0d6SHong Zhang } 1330a42ab0d6SHong Zhang ilen_row++; 1331a42ab0d6SHong Zhang } 1332a42ab0d6SHong Zhang } 1333a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1334a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1335a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1336a42ab0d6SHong Zhang if (tcol) { 1337a42ab0d6SHong Zhang #else 1338a42ab0d6SHong Zhang if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1339a42ab0d6SHong Zhang #endif 1340a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1341a42ab0d6SHong Zhang if (!ijonly) { 1342a42ab0d6SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1343a42ab0d6SHong Zhang mat_a += bs2; 1344a42ab0d6SHong Zhang } 1345a42ab0d6SHong Zhang ilen_row++; 1346a42ab0d6SHong Zhang } 1347a42ab0d6SHong Zhang } 1348a42ab0d6SHong Zhang } else { /* allcolumns */ 1349a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1350a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1351a42ab0d6SHong Zhang *mat_j++ = ctmp; 1352a42ab0d6SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1353a42ab0d6SHong Zhang mat_a += bs2; 1354a42ab0d6SHong Zhang ilen_row++; 1355a42ab0d6SHong Zhang } else break; 1356a42ab0d6SHong Zhang } 1357a42ab0d6SHong Zhang imark = l; 1358a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1359a42ab0d6SHong Zhang *mat_j++ = cstart+cworkA[l]; 1360a42ab0d6SHong Zhang if (!ijonly) { 1361a42ab0d6SHong Zhang ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1362a42ab0d6SHong Zhang mat_a += bs2; 1363a42ab0d6SHong Zhang } 1364a42ab0d6SHong Zhang ilen_row++; 1365a42ab0d6SHong Zhang } 1366a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1367a42ab0d6SHong Zhang *mat_j++ = bmap[cworkB[l]]; 1368a42ab0d6SHong Zhang if (!ijonly) { 1369a42ab0d6SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1370a42ab0d6SHong Zhang mat_a += bs2; 1371a42ab0d6SHong Zhang } 1372a42ab0d6SHong Zhang ilen_row++; 1373a42ab0d6SHong Zhang } 1374a42ab0d6SHong Zhang } 137517df9f7cSHong Zhang imat_ilen[row] = ilen_row; 137617df9f7cSHong Zhang } 137717df9f7cSHong Zhang } 137817df9f7cSHong Zhang } 137917df9f7cSHong Zhang 138017df9f7cSHong Zhang /* Now assemble the off proc rows */ 13815dd0c08cSHong Zhang if (!ijonly) { 138217df9f7cSHong Zhang ierr = MPI_Waitall(nrqs,r_waits4,r_status4);CHKERRQ(ierr); 13835dd0c08cSHong Zhang } 138417df9f7cSHong Zhang for (tmp2=0; tmp2<nrqs; tmp2++) { 138517df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 138617df9f7cSHong Zhang jmax = sbuf1_i[0]; 138717df9f7cSHong Zhang ct1 = 2*jmax + 1; 138817df9f7cSHong Zhang ct2 = 0; 138917df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 139017df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 13915dd0c08cSHong Zhang if (!ijonly) rbuf4_i = rbuf4[tmp2]; 139217df9f7cSHong Zhang for (j=1; j<=jmax; j++) { 139317df9f7cSHong Zhang is_no = sbuf1_i[2*j-1]; 139417df9f7cSHong Zhang rmap_i = rmap[is_no]; 139517df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 139617df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[is_no]->data; 139717df9f7cSHong Zhang imat_ilen = subc->ilen; 139817df9f7cSHong Zhang imat_j = subc->j; 139917df9f7cSHong Zhang imat_i = subc->i; 14005dd0c08cSHong Zhang if (!ijonly) imat_a = subc->a; 140117df9f7cSHong Zhang max1 = sbuf1_i[2*j]; 14025dd0c08cSHong Zhang for (k=0; k<max1; k++,ct1++) { /* for each recved block row */ 140317df9f7cSHong Zhang row = sbuf1_i[ct1]; 140417df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 140517df9f7cSHong Zhang ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 140617df9f7cSHong Zhang row--; 14075dd0c08cSHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 140817df9f7cSHong Zhang #else 140917df9f7cSHong Zhang row = rmap_i[row]; 141017df9f7cSHong Zhang #endif 141117df9f7cSHong Zhang ilen = imat_ilen[row]; 141217df9f7cSHong Zhang mat_i = imat_i[row]; 14135dd0c08cSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 141417df9f7cSHong Zhang mat_j = imat_j + mat_i; 141517df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 141617df9f7cSHong Zhang if (!allcolumns[is_no]) { 141717df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 141817df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 141917df9f7cSHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 142017df9f7cSHong Zhang #else 142117df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 142217df9f7cSHong Zhang #endif 142317df9f7cSHong Zhang if (tcol) { 142417df9f7cSHong Zhang *mat_j++ = tcol - 1; 14255dd0c08cSHong Zhang if (!ijonly) { 14265dd0c08cSHong Zhang ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 14275dd0c08cSHong Zhang mat_a += bs2; 14285dd0c08cSHong Zhang } 14295dd0c08cSHong Zhang //*mat_a++ = rbuf4_i[ct2]; 143017df9f7cSHong Zhang ilen++; 143117df9f7cSHong Zhang } 143217df9f7cSHong Zhang } 143317df9f7cSHong Zhang } else { /* allcolumns */ 143417df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 143517df9f7cSHong Zhang *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 14365dd0c08cSHong Zhang //*mat_a++ = rbuf4_i[ct2]; 14375dd0c08cSHong Zhang if (!ijonly) { 14385dd0c08cSHong Zhang ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 14395dd0c08cSHong Zhang mat_a += bs2; 14405dd0c08cSHong Zhang } 144117df9f7cSHong Zhang ilen++; 144217df9f7cSHong Zhang } 144317df9f7cSHong Zhang } 144417df9f7cSHong Zhang imat_ilen[row] = ilen; 144517df9f7cSHong Zhang } 144617df9f7cSHong Zhang } 144717df9f7cSHong Zhang } 144817df9f7cSHong Zhang 14495dd0c08cSHong Zhang if (!iscsorted) { /* sort column indices of the rows -- not done yet!!! */ 145017df9f7cSHong Zhang for (i=0; i<ismax; i++) { 145117df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 145217df9f7cSHong Zhang imat_j = subc->j; 145317df9f7cSHong Zhang imat_i = subc->i; 145417df9f7cSHong Zhang imat_a = subc->a; 145517df9f7cSHong Zhang imat_ilen = subc->ilen; 145617df9f7cSHong Zhang 145717df9f7cSHong Zhang if (allcolumns[i]) continue; 145817df9f7cSHong Zhang jmax = nrow[i]; 145917df9f7cSHong Zhang for (j=0; j<jmax; j++) { 146017df9f7cSHong Zhang PetscInt ilen; 146117df9f7cSHong Zhang 146217df9f7cSHong Zhang mat_i = imat_i[j]; 14635dd0c08cSHong Zhang if (!ijonly) mat_a = imat_a + mat_i; 146417df9f7cSHong Zhang mat_j = imat_j + mat_i; 146517df9f7cSHong Zhang ilen = imat_ilen[j]; 146617df9f7cSHong Zhang ierr = PetscSortIntWithScalarArray(ilen,mat_j,mat_a);CHKERRQ(ierr); 146717df9f7cSHong Zhang } 146817df9f7cSHong Zhang } 146917df9f7cSHong Zhang } 147017df9f7cSHong Zhang 147117df9f7cSHong Zhang ierr = PetscFree(r_status4);CHKERRQ(ierr); 147217df9f7cSHong Zhang ierr = PetscFree(r_waits4);CHKERRQ(ierr); 147317df9f7cSHong Zhang if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 147417df9f7cSHong Zhang ierr = PetscFree(s_waits4);CHKERRQ(ierr); 147517df9f7cSHong Zhang ierr = PetscFree(s_status4);CHKERRQ(ierr); 147617df9f7cSHong Zhang 147717df9f7cSHong Zhang /* Restore the indices */ 147817df9f7cSHong Zhang for (i=0; i<ismax; i++) { 147917df9f7cSHong Zhang ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 148017df9f7cSHong Zhang if (!allcolumns[i]) { 148117df9f7cSHong Zhang ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 148217df9f7cSHong Zhang } 148317df9f7cSHong Zhang } 148417df9f7cSHong Zhang 148517df9f7cSHong Zhang for (i=0; i<ismax; i++) { 148617df9f7cSHong Zhang ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 148717df9f7cSHong Zhang ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 148817df9f7cSHong Zhang } 148917df9f7cSHong Zhang 149017df9f7cSHong Zhang /* Destroy allocated memory */ 149117df9f7cSHong Zhang if (!ismax) { 149217df9f7cSHong Zhang ierr = PetscFree(pa);CHKERRQ(ierr); 149317df9f7cSHong Zhang 149417df9f7cSHong Zhang ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 149517df9f7cSHong Zhang for (i=0; i<nrqr; ++i) { 149617df9f7cSHong Zhang ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 149717df9f7cSHong Zhang } 149817df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 149917df9f7cSHong Zhang ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 150017df9f7cSHong Zhang } 150117df9f7cSHong Zhang 150217df9f7cSHong Zhang ierr = PetscFree3(sbuf2,req_size,req_source1);CHKERRQ(ierr); 150317df9f7cSHong Zhang ierr = PetscFree3(req_source2,rbuf2,rbuf3);CHKERRQ(ierr); 150417df9f7cSHong Zhang } 150517df9f7cSHong Zhang 150617df9f7cSHong Zhang ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 150717df9f7cSHong Zhang ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 150817df9f7cSHong Zhang ierr = PetscFree5(irow,icol,nrow,ncol,issorted);CHKERRQ(ierr); 150917df9f7cSHong Zhang 151017df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 151117df9f7cSHong Zhang ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 151217df9f7cSHong Zhang } 151317df9f7cSHong Zhang ierr = PetscFree(rbuf4);CHKERRQ(ierr); 151417df9f7cSHong Zhang 151517df9f7cSHong Zhang ierr = PetscFree5(smats,row2proc,cmap,rmap,allcolumns);CHKERRQ(ierr); 151617df9f7cSHong Zhang PetscFunctionReturn(0); 151717df9f7cSHong Zhang } 1518a42ab0d6SHong Zhang 151917df9f7cSHong Zhang //------------ endof new ------------- 152017df9f7cSHong Zhang 15214da72fa9SHong Zhang PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 1522a2feddc5SSatish Balay { 1523df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 1524cf4f527aSSatish Balay Mat A = c->A; 1525df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 15265d0c19d7SBarry Smith const PetscInt **irow,**icol,*irow_i; 15275d0c19d7SBarry Smith PetscInt *nrow,*ncol,*w3,*w4,start; 15286849ba73SBarry Smith PetscErrorCode ierr; 15299405f653SBarry Smith PetscMPIInt size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc; 15309405f653SBarry Smith PetscInt **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row; 1531b24ad042SBarry Smith PetscInt nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol; 1532b24ad042SBarry Smith PetscInt **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 15335d0c19d7SBarry Smith PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 1534052f0c41SBarry Smith PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 1535d0f46423SBarry Smith PetscInt bs =C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 1536899cda47SBarry Smith PetscInt cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 1537899cda47SBarry Smith PetscInt *bmap = c->garray,ctmp,rstart=c->rstartbs; 15387a868f3eSHong Zhang MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3; 15397a868f3eSHong Zhang MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3; 1540d5b485abSSatish Balay MPI_Comm comm; 1541ace3abfcSBarry Smith PetscBool flag; 1542b24ad042SBarry Smith PetscMPIInt *onodes1,*olengths1; 15437a868f3eSHong Zhang PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 154426fbe8dcSKarl Rupp 15457a868f3eSHong Zhang /* variables below are used for the matrix numerical values - case of !ijonly */ 15467a868f3eSHong Zhang MPI_Request *r_waits4,*s_waits4; 15477a868f3eSHong Zhang MPI_Status *r_status4,*s_status4; 15480298fd71SBarry Smith MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a = NULL,*sbuf_aa_i,*vworkA = NULL,*vworkB = NULL; 15497a868f3eSHong Zhang MatScalar *a_a=a->a,*b_a=b->a; 1550c7dd2924SSatish Balay 155180d608c0SSatish Balay #if defined(PETSC_USE_CTABLE) 1552b24ad042SBarry Smith PetscInt tt; 15530298fd71SBarry Smith PetscTable *rmap,*cmap,rmap_i,cmap_i=NULL; 1554233c2ff6SSatish Balay #else 15550298fd71SBarry Smith PetscInt **cmap,*cmap_i=NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 1556233c2ff6SSatish Balay #endif 1557d5b485abSSatish Balay 15583a40ed3dSBarry Smith PetscFunctionBegin; 1559ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 15607adad957SLisandro Dalcin tag0 = ((PetscObject)C)->tag; 1561d5b485abSSatish Balay size = c->size; 1562d5b485abSSatish Balay rank = c->rank; 15639e686edcSHong Zhang //if (!rank) printf("MatGetSubMatrices_MPIBAIJ scall %d, bs %d\n",scall,bs); 1564d5b485abSSatish Balay 156534e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 156634e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 156734e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 156834e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 156934e6ae18SSatish Balay 1570052f0c41SBarry Smith #if defined(PETSC_USE_CTABLE) 1571dcca6d9dSJed Brown ierr = PetscMalloc4(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol);CHKERRQ(ierr); 1572052f0c41SBarry Smith #else 1573dcca6d9dSJed Brown ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,Mbs+1,&rtable);CHKERRQ(ierr); 1574d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 1575d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 157682a7e548SBarry Smith jmax = C->rmap->range[i+1]/bs; 157726fbe8dcSKarl Rupp for (; j<jmax; j++) rtable[j] = i; 1578d5b485abSSatish Balay } 1579233c2ff6SSatish Balay #endif 1580233c2ff6SSatish Balay 1581233c2ff6SSatish Balay for (i=0; i<ismax; i++) { 15824da72fa9SHong Zhang if (allrows[i]) { 15830298fd71SBarry Smith irow[i] = NULL; 15844da72fa9SHong Zhang nrow[i] = C->rmap->N/bs; 15854da72fa9SHong Zhang } else { 1586233c2ff6SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 1587233c2ff6SSatish Balay ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 15884da72fa9SHong Zhang } 15894da72fa9SHong Zhang 1590307b7a18SHong Zhang if (allcolumns[i]) { 15910298fd71SBarry Smith icol[i] = NULL; 1592307b7a18SHong Zhang ncol[i] = C->cmap->N/bs; 1593307b7a18SHong Zhang } else { 1594307b7a18SHong Zhang ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 1595233c2ff6SSatish Balay ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 1596233c2ff6SSatish Balay } 1597307b7a18SHong Zhang } 1598d5b485abSSatish Balay 1599d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg,and buffer space 1600d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 1601884e879aSBarry Smith ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 1602d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1603b24ad042SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 1604d5b485abSSatish Balay jmax = nrow[i]; 1605d5b485abSSatish Balay irow_i = irow[i]; 1606d5b485abSSatish Balay for (j=0; j<jmax; j++) { 160726fbe8dcSKarl Rupp if (allrows[i]) row = j; 160826fbe8dcSKarl Rupp else row = irow_i[j]; 160926fbe8dcSKarl Rupp 1610233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1611899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1612233c2ff6SSatish Balay #else 1613d5b485abSSatish Balay proc = rtable[row]; 1614233c2ff6SSatish Balay #endif 1615d5b485abSSatish Balay w4[proc]++; 1616d5b485abSSatish Balay } 1617d5b485abSSatish Balay for (j=0; j<size; j++) { 1618d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 1619d5b485abSSatish Balay } 1620d5b485abSSatish Balay } 1621d5b485abSSatish Balay 1622d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 1623e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 1624d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 1625d5b485abSSatish Balay w3[rank] = 0; 1626d5b485abSSatish Balay for (i=0; i<size; i++) { 1627d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 1628d5b485abSSatish Balay } 1629854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 1630d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 1631d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 1632d5b485abSSatish Balay } 1633d5b485abSSatish Balay 1634d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 1635d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 1636d5b485abSSatish Balay j = pa[i]; 1637d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 1638d5b485abSSatish Balay msz += w1[j]; 1639d5b485abSSatish Balay } 1640d5b485abSSatish Balay 1641c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 1642a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 1643c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 1644d5b485abSSatish Balay 1645c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 1646c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 1647c7dd2924SSatish Balay 1648c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 1649c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 1650d5b485abSSatish Balay 1651d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 1652dcca6d9dSJed Brown ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 165323bdbc58SBarry Smith ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 165423bdbc58SBarry Smith ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 1655d5b485abSSatish Balay { 1656b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 1657d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 1658d5b485abSSatish Balay j = pa[i]; 1659d5b485abSSatish Balay iptr += ict; 1660d5b485abSSatish Balay sbuf1[j] = iptr; 1661d5b485abSSatish Balay ict = w1[j]; 1662d5b485abSSatish Balay } 1663d5b485abSSatish Balay } 1664d5b485abSSatish Balay 1665d5b485abSSatish Balay /* Form the outgoing messages */ 1666d5b485abSSatish Balay /* Initialise the header space */ 1667d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 1668d5b485abSSatish Balay j = pa[i]; 1669d5b485abSSatish Balay sbuf1[j][0] = 0; 1670b24ad042SBarry Smith ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 1671d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 1672d5b485abSSatish Balay } 1673d5b485abSSatish Balay 1674d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 1675d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1676b24ad042SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 1677d5b485abSSatish Balay irow_i = irow[i]; 1678d5b485abSSatish Balay jmax = nrow[i]; 1679d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 168026fbe8dcSKarl Rupp if (allrows[i]) row = j; 168126fbe8dcSKarl Rupp else row = irow_i[j]; 168226fbe8dcSKarl Rupp 1683233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1684899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1685233c2ff6SSatish Balay #else 1686d5b485abSSatish Balay proc = rtable[row]; 1687233c2ff6SSatish Balay #endif 1688d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 1689d5b485abSSatish Balay ctr[proc]++; 1690d5b485abSSatish Balay *ptr[proc] = row; 1691d5b485abSSatish Balay ptr[proc]++; 1692d5b485abSSatish Balay } 1693d5b485abSSatish Balay } 1694d5b485abSSatish Balay /* Update the headers for the current IS */ 1695d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 169606ef35abSSatish Balay if ((ctr_j = ctr[j])) { 1697d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 1698d5b485abSSatish Balay k = ++sbuf1_j[0]; 1699d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 1700d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 1701d5b485abSSatish Balay } 1702d5b485abSSatish Balay } 1703d5b485abSSatish Balay } 1704d5b485abSSatish Balay 1705d5b485abSSatish Balay /* Now post the sends */ 1706854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 1707d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1708d5b485abSSatish Balay j = pa[i]; 1709b24ad042SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 1710d5b485abSSatish Balay } 1711d5b485abSSatish Balay 1712d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 1713854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 1714854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&rbuf2);CHKERRQ(ierr); 1715d5b485abSSatish Balay rbuf2[0] = tmp + msz; 1716d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 1717d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 1718d5b485abSSatish Balay } 1719d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1720d5b485abSSatish Balay j = pa[i]; 1721b24ad042SBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 1722d5b485abSSatish Balay } 1723d5b485abSSatish Balay 1724d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 1725d5b485abSSatish Balay 1726d5b485abSSatish Balay /* Receive messages*/ 1727854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 1728854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 1729dcca6d9dSJed Brown ierr = PetscMalloc3(nrqr+1,&sbuf2,nrqr,&req_size,nrqr,&req_source);CHKERRQ(ierr); 1730d5b485abSSatish Balay { 1731df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 1732b24ad042SBarry Smith PetscInt *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 1733d5b485abSSatish Balay 1734d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1735999d9058SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 173626fbe8dcSKarl Rupp 1737999d9058SBarry Smith req_size[idex] = 0; 1738999d9058SBarry Smith rbuf1_i = rbuf1[idex]; 1739d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 1740b24ad042SBarry Smith ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 1741785e854fSJed Brown ierr = PetscMalloc1(end,&sbuf2[idex]);CHKERRQ(ierr); 1742999d9058SBarry Smith sbuf2_i = sbuf2[idex]; 1743d5b485abSSatish Balay for (j=start; j<end; j++) { 1744d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 1745d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 1746d5b485abSSatish Balay sbuf2_i[j] = ncols; 1747999d9058SBarry Smith req_size[idex] += ncols; 1748d5b485abSSatish Balay } 1749999d9058SBarry Smith req_source[idex] = r_status1[i].MPI_SOURCE; 1750d5b485abSSatish Balay /* form the header */ 1751999d9058SBarry Smith sbuf2_i[0] = req_size[idex]; 175226fbe8dcSKarl Rupp for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 1753b24ad042SBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 1754d5b485abSSatish Balay } 1755d5b485abSSatish Balay } 1756606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 1757606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 1758d5b485abSSatish Balay 1759d5b485abSSatish Balay /* recv buffer sizes */ 1760d5b485abSSatish Balay /* Receive messages*/ 1761854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&rbuf3);CHKERRQ(ierr); 1762854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 1763854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 17647a868f3eSHong Zhang if (!ijonly) { 1765854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 1766854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 17677a868f3eSHong Zhang } 1768d5b485abSSatish Balay 1769d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1770999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 1771785e854fSJed Brown ierr = PetscMalloc1(rbuf2[idex][0],&rbuf3[idex]);CHKERRQ(ierr); 1772b24ad042SBarry Smith ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 17737a868f3eSHong Zhang if (!ijonly) { 1774785e854fSJed Brown ierr = PetscMalloc1(rbuf2[idex][0]*bs2,&rbuf4[idex]);CHKERRQ(ierr); 1775b24ad042SBarry Smith ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 1776d5b485abSSatish Balay } 17777a868f3eSHong Zhang } 1778606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 1779606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 1780d5b485abSSatish Balay 1781d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 1782854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 1783854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 1784d5b485abSSatish Balay 17850c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 17860c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 1787606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 1788606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 1789606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 1790606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 1791d5b485abSSatish Balay 1792d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 1793854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 1794d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1795854ce69bSBarry Smith ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 1796d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 1797d5b485abSSatish Balay 1798854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 1799d5b485abSSatish Balay { 1800d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1801d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1802d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 1803d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 1804d5b485abSSatish Balay ct2 = 0; 1805d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1806d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 1807d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1808e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 180905b2c859SBarry Smith nzA = a_i[row+1] - a_i[row]; 181005b2c859SBarry Smith nzB = b_i[row+1] - b_i[row]; 1811d5b485abSSatish Balay ncols = nzA + nzB; 181205b2c859SBarry Smith cworkA = a_j + a_i[row]; 181305b2c859SBarry Smith cworkB = b_j + b_i[row]; 1814d5b485abSSatish Balay 1815d5b485abSSatish Balay /* load the column indices for this row into cols*/ 1816d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 1817d5b485abSSatish Balay for (l=0; l<nzB; l++) { 1818d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 1819d5b485abSSatish Balay else break; 1820d5b485abSSatish Balay } 1821d5b485abSSatish Balay imark = l; 1822d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 1823d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 1824d5b485abSSatish Balay ct2 += ncols; 1825d5b485abSSatish Balay } 1826d5b485abSSatish Balay } 1827b24ad042SBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 1828d5b485abSSatish Balay } 1829d5b485abSSatish Balay } 1830854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_status3);CHKERRQ(ierr); 1831854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_status3);CHKERRQ(ierr); 1832d5b485abSSatish Balay 1833d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 18347a868f3eSHong Zhang if (!ijonly) { 1835854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 1836d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1837785e854fSJed Brown ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 1838a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1839d5b485abSSatish Balay 1840854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 1841d5b485abSSatish Balay { 1842d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1843d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1844d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 1845d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 1846d5b485abSSatish Balay ct2 = 0; 1847d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1848d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 1849d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1850e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 185105b2c859SBarry Smith nzA = a_i[row+1] - a_i[row]; 185205b2c859SBarry Smith nzB = b_i[row+1] - b_i[row]; 1853d5b485abSSatish Balay ncols = nzA + nzB; 185405b2c859SBarry Smith cworkB = b_j + b_i[row]; 185505b2c859SBarry Smith vworkA = a_a + a_i[row]*bs2; 185605b2c859SBarry Smith vworkB = b_a + b_i[row]*bs2; 1857d5b485abSSatish Balay 1858d5b485abSSatish Balay /* load the column values for this row into vals*/ 18595b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 1860d5b485abSSatish Balay for (l=0; l<nzB; l++) { 18613a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 18623eda8832SBarry Smith ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 186326fbe8dcSKarl Rupp } else break; 1864d5b485abSSatish Balay } 1865d5b485abSSatish Balay imark = l; 18663a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 18673eda8832SBarry Smith ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 18683a40ed3dSBarry Smith } 18693a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 18703eda8832SBarry Smith ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 18713a40ed3dSBarry Smith } 1872d5b485abSSatish Balay ct2 += ncols; 1873d5b485abSSatish Balay } 1874d5b485abSSatish Balay } 18753eda8832SBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 1876d5b485abSSatish Balay } 1877d5b485abSSatish Balay } 1878854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 1879854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 18807a868f3eSHong Zhang } 1881533163c2SBarry Smith ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 1882606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1883d5b485abSSatish Balay 1884d5b485abSSatish Balay /* Form the matrix */ 1885307b7a18SHong Zhang /* create col map: global col of C -> local col of submatrices */ 1886d5b485abSSatish Balay { 18875d0c19d7SBarry Smith const PetscInt *icol_i; 1888233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1889854ce69bSBarry Smith ierr = PetscMalloc1(1+ismax,&cmap);CHKERRQ(ierr); 1890ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 1891307b7a18SHong Zhang if (!allcolumns[i]) { 18928fa52d88SHong Zhang ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 1893307b7a18SHong Zhang jmax = ncol[i]; 1894307b7a18SHong Zhang icol_i = icol[i]; 18958fa52d88SHong Zhang cmap_i = cmap[i]; 1896307b7a18SHong Zhang for (j=0; j<jmax; j++) { 18973861aac3SJed Brown ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1898307b7a18SHong Zhang } 1899307b7a18SHong Zhang } else { 19000298fd71SBarry Smith cmap[i] = NULL; 1901307b7a18SHong Zhang } 1902233c2ff6SSatish Balay } 1903233c2ff6SSatish Balay #else 1904785e854fSJed Brown ierr = PetscMalloc1(ismax,&cmap);CHKERRQ(ierr); 1905d5b485abSSatish Balay for (i=0; i<ismax; i++) { 19068fa52d88SHong Zhang if (!allcolumns[i]) { 1907884e879aSBarry Smith ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 1908d5b485abSSatish Balay jmax = ncol[i]; 1909d5b485abSSatish Balay icol_i = icol[i]; 1910d5b485abSSatish Balay cmap_i = cmap[i]; 1911d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1912d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1913d5b485abSSatish Balay } 19148fa52d88SHong Zhang } else { /* allcolumns[i] */ 19150298fd71SBarry Smith cmap[i] = NULL; 19168fa52d88SHong Zhang } 1917d5b485abSSatish Balay } 1918307b7a18SHong Zhang #endif 1919d5b485abSSatish Balay } 1920d5b485abSSatish Balay 1921d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 192226fbe8dcSKarl Rupp for (i=0,j=0; i<ismax; i++) j += nrow[i]; 1923884e879aSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr); 1924b24ad042SBarry Smith lens[0] = (PetscInt*)(lens + ismax); 1925b24ad042SBarry Smith ierr = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr); 192626fbe8dcSKarl Rupp for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 1927d5b485abSSatish Balay 1928d5b485abSSatish Balay /* Update lens from local data */ 1929d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1930d5b485abSSatish Balay jmax = nrow[i]; 19318fa52d88SHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 1932d5b485abSSatish Balay irow_i = irow[i]; 1933d5b485abSSatish Balay lens_i = lens[i]; 1934d5b485abSSatish Balay for (j=0; j<jmax; j++) { 193526fbe8dcSKarl Rupp if (allrows[i]) row = j; 193626fbe8dcSKarl Rupp else row = irow_i[j]; 193726fbe8dcSKarl Rupp 1938233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1939899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1940233c2ff6SSatish Balay #else 1941d5b485abSSatish Balay proc = rtable[row]; 1942233c2ff6SSatish Balay #endif 1943d5b485abSSatish Balay if (proc == rank) { 194436f4e84dSSatish Balay /* Get indices from matA and then from matB */ 194536f4e84dSSatish Balay row = row - rstart; 194605b2c859SBarry Smith nzA = a_i[row+1] - a_i[row]; 194705b2c859SBarry Smith nzB = b_i[row+1] - b_i[row]; 194805b2c859SBarry Smith cworkA = a_j + a_i[row]; 194905b2c859SBarry Smith cworkB = b_j + b_i[row]; 1950307b7a18SHong Zhang if (!allcolumns[i]) { 1951233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1952233c2ff6SSatish Balay for (k=0; k<nzA; k++) { 19538fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 195426fbe8dcSKarl Rupp if (tt) lens_i[j]++; 1955233c2ff6SSatish Balay } 1956233c2ff6SSatish Balay for (k=0; k<nzB; k++) { 19578fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 195826fbe8dcSKarl Rupp if (tt) lens_i[j]++; 1959233c2ff6SSatish Balay } 1960307b7a18SHong Zhang 1961233c2ff6SSatish Balay #else 1962ca161407SBarry Smith for (k=0; k<nzA; k++) { 196326fbe8dcSKarl Rupp if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1964ca161407SBarry Smith } 1965ca161407SBarry Smith for (k=0; k<nzB; k++) { 196626fbe8dcSKarl Rupp if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1967d5b485abSSatish Balay } 1968233c2ff6SSatish Balay #endif 1969307b7a18SHong Zhang } else { /* allcolumns */ 1970307b7a18SHong Zhang lens_i[j] = nzA + nzB; 1971307b7a18SHong Zhang } 1972d5b485abSSatish Balay } 1973a2feddc5SSatish Balay } 1974ca161407SBarry Smith } 1975233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 1976d5b485abSSatish Balay /* Create row map*/ 1977854ce69bSBarry Smith ierr = PetscMalloc1(1+ismax,&rmap);CHKERRQ(ierr); 1978ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 19798fa52d88SHong Zhang ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1980233c2ff6SSatish Balay } 1981233c2ff6SSatish Balay #else 1982233c2ff6SSatish Balay /* Create row map*/ 1983884e879aSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr); 1984b24ad042SBarry Smith rmap[0] = (PetscInt*)(rmap + ismax); 1985b24ad042SBarry Smith ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr); 198626fbe8dcSKarl Rupp for (i=1; i<ismax; i++) rmap[i] = rmap[i-1] + Mbs; 1987233c2ff6SSatish Balay #endif 1988d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1989d5b485abSSatish Balay irow_i = irow[i]; 1990d5b485abSSatish Balay jmax = nrow[i]; 1991233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 19928fa52d88SHong Zhang rmap_i = rmap[i]; 1993233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 19944da72fa9SHong Zhang if (allrows[i]) { 19953861aac3SJed Brown ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 19964da72fa9SHong Zhang } else { 19973861aac3SJed Brown ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1998233c2ff6SSatish Balay } 19994da72fa9SHong Zhang } 2000233c2ff6SSatish Balay #else 2001233c2ff6SSatish Balay rmap_i = rmap[i]; 2002d5b485abSSatish Balay for (j=0; j<jmax; j++) { 200326fbe8dcSKarl Rupp if (allrows[i]) rmap_i[j] = j; 200426fbe8dcSKarl Rupp else rmap_i[irow_i[j]] = j; 20054da72fa9SHong Zhang } 2006233c2ff6SSatish Balay #endif 2007d5b485abSSatish Balay } 2008d5b485abSSatish Balay 2009d5b485abSSatish Balay /* Update lens from offproc data */ 2010d5b485abSSatish Balay { 2011b24ad042SBarry Smith PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 2012b24ad042SBarry Smith PetscMPIInt ii; 2013d5b485abSSatish Balay 2014d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 2015b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr); 2016b24ad042SBarry Smith idex = pa[ii]; 2017999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 2018d5b485abSSatish Balay jmax = sbuf1_i[0]; 2019d5b485abSSatish Balay ct1 = 2*jmax+1; 2020d5b485abSSatish Balay ct2 = 0; 2021b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 2022b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 2023d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 2024d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 2025d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 2026d5b485abSSatish Balay lens_i = lens[is_no]; 20278fa52d88SHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 2028d5b485abSSatish Balay rmap_i = rmap[is_no]; 2029d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 2030233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 20318fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 2032233c2ff6SSatish Balay row--; 203326fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2034233c2ff6SSatish Balay #else 2035d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 2036233c2ff6SSatish Balay #endif 2037d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 2038d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 2039307b7a18SHong Zhang if (!allcolumns[is_no]) { 2040233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 20418fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 204226fbe8dcSKarl Rupp if (tt) lens_i[row]++; 2043233c2ff6SSatish Balay #else 204426fbe8dcSKarl Rupp if (cmap_i[rbuf3_i[ct2]]) lens_i[row]++; 2045233c2ff6SSatish Balay #endif 2046307b7a18SHong Zhang } else { /* allcolumns */ 2047307b7a18SHong Zhang lens_i[row]++; 2048307b7a18SHong Zhang } 2049d5b485abSSatish Balay } 2050d5b485abSSatish Balay } 2051d5b485abSSatish Balay } 2052d5b485abSSatish Balay } 2053d5b485abSSatish Balay } 2054606d414cSSatish Balay ierr = PetscFree(r_status3);CHKERRQ(ierr); 2055606d414cSSatish Balay ierr = PetscFree(r_waits3);CHKERRQ(ierr); 20560c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 2057606d414cSSatish Balay ierr = PetscFree(s_status3);CHKERRQ(ierr); 2058606d414cSSatish Balay ierr = PetscFree(s_waits3);CHKERRQ(ierr); 2059d5b485abSSatish Balay 2060d5b485abSSatish Balay /* Create the submatrices */ 2061d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 20627a868f3eSHong Zhang if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet"); 2063d5b485abSSatish Balay /* 2064d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 2065d5b485abSSatish Balay */ 2066d5b485abSSatish Balay for (i=0; i<ismax; i++) { 2067df36731bSBarry Smith mat = (Mat_SeqBAIJ*)(submats[i]->data); 2068e7e72b3dSBarry 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"); 2069b24ad042SBarry Smith ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 2070e7e72b3dSBarry Smith if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 2071d5b485abSSatish Balay /* Initial matrix as if empty */ 2072b24ad042SBarry Smith ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr); 207326fbe8dcSKarl Rupp 2074d5f3da31SBarry Smith submats[i]->factortype = C->factortype; 2075d5b485abSSatish Balay } 2076ca161407SBarry Smith } else { 20777a868f3eSHong Zhang PetscInt bs_tmp; 207826fbe8dcSKarl Rupp if (ijonly) bs_tmp = 1; 207926fbe8dcSKarl Rupp else bs_tmp = bs; 2080d5b485abSSatish Balay for (i=0; i<ismax; i++) { 2081f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 20827a868f3eSHong Zhang ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr); 20837adad957SLisandro Dalcin ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 20847a868f3eSHong Zhang ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 20857a868f3eSHong Zhang ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 2086d5b485abSSatish Balay } 2087d5b485abSSatish Balay } 2088d5b485abSSatish Balay 2089d5b485abSSatish Balay /* Assemble the matrices */ 2090d5b485abSSatish Balay /* First assemble the local rows */ 2091d5b485abSSatish Balay { 2092b24ad042SBarry Smith PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 20930298fd71SBarry Smith MatScalar *imat_a = NULL; 2094d5b485abSSatish Balay 2095d5b485abSSatish Balay for (i=0; i<ismax; i++) { 2096df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[i]->data; 2097d5b485abSSatish Balay imat_ilen = mat->ilen; 2098d5b485abSSatish Balay imat_j = mat->j; 2099d5b485abSSatish Balay imat_i = mat->i; 21007a868f3eSHong Zhang if (!ijonly) imat_a = mat->a; 21018fa52d88SHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 2102d5b485abSSatish Balay rmap_i = rmap[i]; 2103d5b485abSSatish Balay irow_i = irow[i]; 2104d5b485abSSatish Balay jmax = nrow[i]; 2105d5b485abSSatish Balay for (j=0; j<jmax; j++) { 210626fbe8dcSKarl Rupp if (allrows[i]) row = j; 210726fbe8dcSKarl Rupp else row = irow_i[j]; 210826fbe8dcSKarl Rupp 2109233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 2110899cda47SBarry Smith ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 2111233c2ff6SSatish Balay #else 2112d5b485abSSatish Balay proc = rtable[row]; 2113233c2ff6SSatish Balay #endif 2114d5b485abSSatish Balay if (proc == rank) { 211536f4e84dSSatish Balay row = row - rstart; 211636f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 211736f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 211836f4e84dSSatish Balay cworkA = a_j + a_i[row]; 211936f4e84dSSatish Balay cworkB = b_j + b_i[row]; 21207a868f3eSHong Zhang if (!ijonly) { 212136f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 212236f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 21237a868f3eSHong Zhang } 2124233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 21258fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 2126233c2ff6SSatish Balay row--; 212726fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2128233c2ff6SSatish Balay #else 212936f4e84dSSatish Balay row = rmap_i[row + rstart]; 2130233c2ff6SSatish Balay #endif 2131df36731bSBarry Smith mat_i = imat_i[row]; 21327a868f3eSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 2133d5b485abSSatish Balay mat_j = imat_j + mat_i; 213436f4e84dSSatish Balay ilen_row = imat_ilen[row]; 213536f4e84dSSatish Balay 213636f4e84dSSatish Balay /* load the column indices for this row into cols*/ 2137307b7a18SHong Zhang if (!allcolumns[i]) { 213836f4e84dSSatish Balay for (l=0; l<nzB; l++) { 213936f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 2140233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 21418fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 2142233c2ff6SSatish Balay if (tcol) { 2143233c2ff6SSatish Balay #else 214436f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 2145233c2ff6SSatish Balay #endif 2146df36731bSBarry Smith *mat_j++ = tcol - 1; 21473eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2148549d3d68SSatish Balay mat_a += bs2; 2149d5b485abSSatish Balay ilen_row++; 2150d5b485abSSatish Balay } 2151ca161407SBarry Smith } else break; 215236f4e84dSSatish Balay } 215336f4e84dSSatish Balay imark = l; 215436f4e84dSSatish Balay for (l=0; l<nzA; l++) { 2155233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 21568fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 2157233c2ff6SSatish Balay if (tcol) { 2158233c2ff6SSatish Balay #else 215936f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 2160233c2ff6SSatish Balay #endif 216136f4e84dSSatish Balay *mat_j++ = tcol - 1; 21627a868f3eSHong Zhang if (!ijonly) { 21633eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2164549d3d68SSatish Balay mat_a += bs2; 21657a868f3eSHong Zhang } 216636f4e84dSSatish Balay ilen_row++; 216736f4e84dSSatish Balay } 216836f4e84dSSatish Balay } 216936f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 2170233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 21718fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 2172233c2ff6SSatish Balay if (tcol) { 2173233c2ff6SSatish Balay #else 217436f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 2175233c2ff6SSatish Balay #endif 217636f4e84dSSatish Balay *mat_j++ = tcol - 1; 21777a868f3eSHong Zhang if (!ijonly) { 21783eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2179549d3d68SSatish Balay mat_a += bs2; 21807a868f3eSHong Zhang } 218136f4e84dSSatish Balay ilen_row++; 218236f4e84dSSatish Balay } 218336f4e84dSSatish Balay } 2184307b7a18SHong Zhang } else { /* allcolumns */ 2185307b7a18SHong Zhang for (l=0; l<nzB; l++) { 2186307b7a18SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 2187307b7a18SHong Zhang *mat_j++ = ctmp; 2188307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2189307b7a18SHong Zhang mat_a += bs2; 2190307b7a18SHong Zhang ilen_row++; 2191307b7a18SHong Zhang } else break; 2192307b7a18SHong Zhang } 2193307b7a18SHong Zhang imark = l; 2194307b7a18SHong Zhang for (l=0; l<nzA; l++) { 2195307b7a18SHong Zhang *mat_j++ = cstart+cworkA[l]; 2196307b7a18SHong Zhang if (!ijonly) { 2197307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2198307b7a18SHong Zhang mat_a += bs2; 2199307b7a18SHong Zhang } 2200307b7a18SHong Zhang ilen_row++; 2201307b7a18SHong Zhang } 2202307b7a18SHong Zhang for (l=imark; l<nzB; l++) { 2203307b7a18SHong Zhang *mat_j++ = bmap[cworkB[l]]; 2204307b7a18SHong Zhang if (!ijonly) { 2205307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2206307b7a18SHong Zhang mat_a += bs2; 2207307b7a18SHong Zhang } 2208307b7a18SHong Zhang ilen_row++; 2209307b7a18SHong Zhang } 2210307b7a18SHong Zhang } 2211d5b485abSSatish Balay imat_ilen[row] = ilen_row; 2212d5b485abSSatish Balay } 2213d5b485abSSatish Balay } 2214d5b485abSSatish Balay } 2215d5b485abSSatish Balay } 2216d5b485abSSatish Balay 2217d5b485abSSatish Balay /* Now assemble the off proc rows*/ 2218d5b485abSSatish Balay { 2219b24ad042SBarry Smith PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 2220b24ad042SBarry Smith PetscInt *imat_j,*imat_i; 22210298fd71SBarry Smith MatScalar *imat_a = NULL,*rbuf4_i = NULL; 2222b24ad042SBarry Smith PetscMPIInt ii; 2223d5b485abSSatish Balay 2224d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 222526fbe8dcSKarl Rupp if (ijonly) ii = tmp2; 222626fbe8dcSKarl Rupp else { 2227b24ad042SBarry Smith ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr); 22287a868f3eSHong Zhang } 2229b24ad042SBarry Smith idex = pa[ii]; 2230999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 2231d5b485abSSatish Balay jmax = sbuf1_i[0]; 2232d5b485abSSatish Balay ct1 = 2*jmax + 1; 2233d5b485abSSatish Balay ct2 = 0; 2234b24ad042SBarry Smith rbuf2_i = rbuf2[ii]; 2235b24ad042SBarry Smith rbuf3_i = rbuf3[ii]; 22367a868f3eSHong Zhang if (!ijonly) rbuf4_i = rbuf4[ii]; 2237d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 2238d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 22398fa52d88SHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 2240d5b485abSSatish Balay rmap_i = rmap[is_no]; 2241df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[is_no]->data; 2242d5b485abSSatish Balay imat_ilen = mat->ilen; 2243d5b485abSSatish Balay imat_j = mat->j; 2244d5b485abSSatish Balay imat_i = mat->i; 22457a868f3eSHong Zhang if (!ijonly) imat_a = mat->a; 2246d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 2247d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 2248d5b485abSSatish Balay row = sbuf1_i[ct1]; 2249233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 22508fa52d88SHong Zhang ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 2251233c2ff6SSatish Balay row--; 225226fbe8dcSKarl Rupp if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2253233c2ff6SSatish Balay #else 2254d5b485abSSatish Balay row = rmap_i[row]; 2255233c2ff6SSatish Balay #endif 2256d5b485abSSatish Balay ilen = imat_ilen[row]; 2257df36731bSBarry Smith mat_i = imat_i[row]; 22587a868f3eSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 2259d5b485abSSatish Balay mat_j = imat_j + mat_i; 2260d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 2261307b7a18SHong Zhang 2262307b7a18SHong Zhang if (!allcolumns[is_no]) { 2263d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 2264233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 22658fa52d88SHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 2266233c2ff6SSatish Balay if (tcol) { 2267233c2ff6SSatish Balay #else 2268d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 2269233c2ff6SSatish Balay #endif 2270df36731bSBarry Smith *mat_j++ = tcol - 1; 22717a868f3eSHong Zhang if (!ijonly) { 22723eda8832SBarry Smith ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2273549d3d68SSatish Balay mat_a += bs2; 22747a868f3eSHong Zhang } 2275d5b485abSSatish Balay ilen++; 2276d5b485abSSatish Balay } 2277d5b485abSSatish Balay } 2278307b7a18SHong Zhang } else { /* allcolumns */ 2279307b7a18SHong Zhang for (l=0; l<max2; l++,ct2++) { 2280307b7a18SHong Zhang *mat_j++ = rbuf3_i[ct2]; 2281307b7a18SHong Zhang if (!ijonly) { 2282307b7a18SHong Zhang ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2283307b7a18SHong Zhang mat_a += bs2; 2284307b7a18SHong Zhang } 2285307b7a18SHong Zhang ilen++; 2286307b7a18SHong Zhang } 2287307b7a18SHong Zhang } 2288d5b485abSSatish Balay imat_ilen[row] = ilen; 2289d5b485abSSatish Balay } 2290d5b485abSSatish Balay } 2291d5b485abSSatish Balay } 2292d5b485abSSatish Balay } 2293cdc6f3adSToby Isaac /* sort the rows */ 2294cdc6f3adSToby Isaac { 2295cdc6f3adSToby Isaac PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 2296cdc6f3adSToby Isaac MatScalar *imat_a = NULL; 2297cdc6f3adSToby Isaac MatScalar *work; 2298cdc6f3adSToby Isaac 2299cdc6f3adSToby Isaac ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 2300cdc6f3adSToby Isaac for (i=0; i<ismax; i++) { 2301cdc6f3adSToby Isaac mat = (Mat_SeqBAIJ*)submats[i]->data; 2302cdc6f3adSToby Isaac imat_ilen = mat->ilen; 2303cdc6f3adSToby Isaac imat_j = mat->j; 2304cdc6f3adSToby Isaac imat_i = mat->i; 2305cdc6f3adSToby Isaac if (!ijonly) imat_a = mat->a; 2306cdc6f3adSToby Isaac if (allcolumns[i]) continue; 2307cdc6f3adSToby Isaac jmax = nrow[i]; 2308cdc6f3adSToby Isaac for (j=0; j<jmax; j++) { 2309cdc6f3adSToby Isaac mat_i = imat_i[j]; 2310cdc6f3adSToby Isaac if (!ijonly) mat_a = imat_a + mat_i*bs2; 2311cdc6f3adSToby Isaac mat_j = imat_j + mat_i; 2312cdc6f3adSToby Isaac ilen_row = imat_ilen[j]; 2313cdc6f3adSToby Isaac if (!ijonly) {ierr = PetscSortIntWithDataArray(ilen_row,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);} 2314cdc6f3adSToby Isaac else {ierr = PetscSortInt(ilen_row,mat_j);CHKERRQ(ierr);} 2315cdc6f3adSToby Isaac } 2316cdc6f3adSToby Isaac } 2317cdc6f3adSToby Isaac ierr = PetscFree(work);CHKERRQ(ierr); 2318cdc6f3adSToby Isaac } 23197a868f3eSHong Zhang if (!ijonly) { 2320606d414cSSatish Balay ierr = PetscFree(r_status4);CHKERRQ(ierr); 2321606d414cSSatish Balay ierr = PetscFree(r_waits4);CHKERRQ(ierr); 23220c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 2323606d414cSSatish Balay ierr = PetscFree(s_waits4);CHKERRQ(ierr); 2324606d414cSSatish Balay ierr = PetscFree(s_status4);CHKERRQ(ierr); 23257a868f3eSHong Zhang } 2326d5b485abSSatish Balay 2327d5b485abSSatish Balay /* Restore the indices */ 2328d5b485abSSatish Balay for (i=0; i<ismax; i++) { 23294da72fa9SHong Zhang if (!allrows[i]) { 2330d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 23314da72fa9SHong Zhang } 2332307b7a18SHong Zhang if (!allcolumns[i]) { 2333d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 2334d5b485abSSatish Balay } 2335307b7a18SHong Zhang } 2336d5b485abSSatish Balay 2337d5b485abSSatish Balay /* Destroy allocated memory */ 233823bdbc58SBarry Smith #if defined(PETSC_USE_CTABLE) 233923bdbc58SBarry Smith ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr); 234023bdbc58SBarry Smith #else 234123bdbc58SBarry Smith ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr); 234223bdbc58SBarry Smith #endif 234323bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 2344606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 2345d5b485abSSatish Balay 234623bdbc58SBarry Smith ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 2347606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 2348606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 2349d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 2350606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 2351d5b485abSSatish Balay } 2352d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 2353606d414cSSatish Balay ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 2354d5b485abSSatish Balay } 235523bdbc58SBarry Smith ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr); 2356606d414cSSatish Balay ierr = PetscFree(rbuf3);CHKERRQ(ierr); 2357606d414cSSatish Balay ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 2358606d414cSSatish Balay ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 23597a868f3eSHong Zhang if (!ijonly) { 23607a868f3eSHong Zhang for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);} 23617a868f3eSHong Zhang ierr = PetscFree(rbuf4);CHKERRQ(ierr); 2362606d414cSSatish Balay ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 2363606d414cSSatish Balay ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 23647a868f3eSHong Zhang } 2365d5b485abSSatish Balay 2366233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 2367ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 23688fa52d88SHong Zhang ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr); 2369233c2ff6SSatish Balay } 2370233c2ff6SSatish Balay #endif 23718fa52d88SHong Zhang ierr = PetscFree(rmap);CHKERRQ(ierr); 23728fa52d88SHong Zhang 23738fa52d88SHong Zhang for (i=0; i<ismax; i++) { 23748fa52d88SHong Zhang if (!allcolumns[i]) { 23758fa52d88SHong Zhang #if defined(PETSC_USE_CTABLE) 23768fa52d88SHong Zhang ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr); 23778fa52d88SHong Zhang #else 23788fa52d88SHong Zhang ierr = PetscFree(cmap[i]);CHKERRQ(ierr); 23798fa52d88SHong Zhang #endif 23808fa52d88SHong Zhang } 23818fa52d88SHong Zhang } 23828fa52d88SHong Zhang ierr = PetscFree(cmap);CHKERRQ(ierr); 2383606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 2384d5b485abSSatish Balay 2385d5b485abSSatish Balay for (i=0; i<ismax; i++) { 238636f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 238736f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2388d5b485abSSatish Balay } 23897a868f3eSHong Zhang 23907a868f3eSHong Zhang c->ijonly = PETSC_FALSE; /* set back to the default */ 23913a40ed3dSBarry Smith PetscFunctionReturn(0); 2392d5b485abSSatish Balay } 2393ca161407SBarry Smith 2394