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; 69749130a8SStefano Zampini MPI_Comm comm,*iscomms; 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 83f489ac74SBarry Smith ierr = PetscMalloc2(imax+1,(PetscInt***)&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++) { 94580bdb30SBarry Smith ierr = PetscArrayzero(w4,size);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); 137580bdb30SBarry Smith ierr = PetscArrayzero(outdat,size);CHKERRQ(ierr); 138580bdb30SBarry Smith ierr = PetscArrayzero(ptr,size);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; 154580bdb30SBarry Smith ierr = PetscArrayzero(outdat[j]+1,2*w3[j]);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++) { 174580bdb30SBarry Smith ierr = PetscArrayzero(ctr,size);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 } 2152cd485c2SBarry Smith ierr = PetscFree2(*(PetscInt***)&idx,n);CHKERRQ(ierr); 216d5b485abSSatish Balay 217749130a8SStefano Zampini ierr = PetscMalloc1(imax,&iscomms);CHKERRQ(ierr); 218d5b485abSSatish Balay for (i=0; i<imax; ++i) { 219749130a8SStefano Zampini ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]),&iscomms[i],NULL);CHKERRQ(ierr); 2206bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 221d5b485abSSatish Balay } 222d5b485abSSatish Balay 223d5b485abSSatish Balay /* Do Local work*/ 224df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 225d5b485abSSatish Balay 226d5b485abSSatish Balay /* Receive messages*/ 227854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&recv_status);CHKERRQ(ierr); 2280c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);} 229d5b485abSSatish Balay 230854ce69bSBarry Smith ierr = PetscMalloc1(nrqs+1,&s_status);CHKERRQ(ierr); 2310c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 232d5b485abSSatish Balay 233d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 23423bdbc58SBarry Smith ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr); 23523bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 236d5b485abSSatish Balay 237854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&xdata);CHKERRQ(ierr); 238854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&isz1);CHKERRQ(ierr); 239df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 240c05d87d6SBarry Smith ierr = PetscFree(rbuf[0]);CHKERRQ(ierr); 241606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 242d5b485abSSatish Balay 243d5b485abSSatish Balay /* Send the data back*/ 244d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 245d5b485abSSatish Balay { 246b24ad042SBarry Smith PetscMPIInt *rw1; 247d5b485abSSatish Balay 248884e879aSBarry Smith ierr = PetscCalloc1(size,&rw1);CHKERRQ(ierr); 249c7dd2924SSatish Balay 250d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 251d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 252e32f2f54SBarry Smith if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 253d5b485abSSatish Balay rw1[proc] = isz1[i]; 254d5b485abSSatish Balay } 255d5b485abSSatish Balay 256c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 257c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 258d5b485abSSatish Balay 259c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 260c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 26103512dc5SSatish Balay ierr = PetscFree(rw1);CHKERRQ(ierr); 262c7dd2924SSatish Balay } 263c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 264c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 265d5b485abSSatish Balay 266d5b485abSSatish Balay /* Now post the sends */ 267854ce69bSBarry Smith ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 268d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 269d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 270b24ad042SBarry Smith ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 271d5b485abSSatish Balay } 272d5b485abSSatish Balay 273d5b485abSSatish Balay /* receive work done on other processors*/ 274d5b485abSSatish Balay { 275b24ad042SBarry Smith PetscMPIInt idex; 276b24ad042SBarry Smith PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 277f1af5d2fSBarry Smith PetscBT table_i; 278d5b485abSSatish Balay MPI_Status *status2; 279d5b485abSSatish Balay 280854ce69bSBarry Smith ierr = PetscMalloc1(PetscMax(nrqr,nrqs)+1,&status2);CHKERRQ(ierr); 281d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 282999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 283d5b485abSSatish Balay /* Process the message*/ 284999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 285d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 286999d9058SBarry Smith jmax = rbuf2[idex][0]; 287d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 288d5b485abSSatish Balay max = rbuf2_i[2*j]; 289d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 290d5b485abSSatish Balay isz_i = isz[is_no]; 291d5b485abSSatish Balay data_i = data[is_no]; 292d5b485abSSatish Balay table_i = table[is_no]; 293d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 294d5b485abSSatish Balay row = rbuf2_i[ct1]; 29526fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 296d5b485abSSatish Balay } 297d5b485abSSatish Balay isz[is_no] = isz_i; 298d5b485abSSatish Balay } 299d5b485abSSatish Balay } 3000c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);} 301606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 302d5b485abSSatish Balay } 303d5b485abSSatish Balay 304d5b485abSSatish Balay for (i=0; i<imax; ++i) { 305749130a8SStefano Zampini ierr = ISCreateGeneral(iscomms[i],isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 306749130a8SStefano Zampini ierr = PetscCommDestroy(&iscomms[i]);CHKERRQ(ierr); 307d5b485abSSatish Balay } 308d5b485abSSatish Balay 309749130a8SStefano Zampini ierr = PetscFree(iscomms);CHKERRQ(ierr); 310c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 311c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 312c7dd2924SSatish Balay 313606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 314c05d87d6SBarry Smith ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr); 315606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 316606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 317606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 318606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 319606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 3208f9f447aSHong Zhang ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr); 321606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 322606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 323606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 324606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 325606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 3263a40ed3dSBarry Smith PetscFunctionReturn(0); 327d5b485abSSatish Balay } 328d5b485abSSatish Balay 329d5b485abSSatish Balay /* 330df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 331d5b485abSSatish Balay the work on the local processor. 332d5b485abSSatish Balay 333d5b485abSSatish Balay Inputs: 334df36731bSBarry Smith C - MAT_MPIBAIJ; 335d5b485abSSatish Balay imax - total no of index sets processed at a time; 336df36731bSBarry Smith table - an array of char - size = Mbs bits. 337d5b485abSSatish Balay 338d5b485abSSatish Balay Output: 3390e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 340d5b485abSSatish Balay to each index set; 341d5b485abSSatish Balay data - pointer to the solutions 342d5b485abSSatish Balay */ 343b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 344d5b485abSSatish Balay { 345df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 346d5b485abSSatish Balay Mat A = c->A,B = c->B; 347df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 348b24ad042SBarry Smith PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 349b24ad042SBarry Smith PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 350f1af5d2fSBarry Smith PetscBT table_i; 351d5b485abSSatish Balay 3523a40ed3dSBarry Smith PetscFunctionBegin; 353899cda47SBarry Smith rstart = c->rstartbs; 354899cda47SBarry Smith cstart = c->cstartbs; 355d5b485abSSatish Balay ai = a->i; 356df36731bSBarry Smith aj = a->j; 357d5b485abSSatish Balay bi = b->i; 358df36731bSBarry Smith bj = b->j; 359d5b485abSSatish Balay garray = c->garray; 360d5b485abSSatish Balay 361d5b485abSSatish Balay 362d5b485abSSatish Balay for (i=0; i<imax; i++) { 363d5b485abSSatish Balay data_i = data[i]; 364d5b485abSSatish Balay table_i = table[i]; 365d5b485abSSatish Balay isz_i = isz[i]; 366d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 367d5b485abSSatish Balay row = data_i[j] - rstart; 368d5b485abSSatish Balay start = ai[row]; 369d5b485abSSatish Balay end = ai[row+1]; 370d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 371df36731bSBarry Smith val = aj[k] + cstart; 37226fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 373d5b485abSSatish Balay } 374d5b485abSSatish Balay start = bi[row]; 375d5b485abSSatish Balay end = bi[row+1]; 376d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 377df36731bSBarry Smith val = garray[bj[k]]; 37826fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 379d5b485abSSatish Balay } 380d5b485abSSatish Balay } 381d5b485abSSatish Balay isz[i] = isz_i; 382d5b485abSSatish Balay } 3833a40ed3dSBarry Smith PetscFunctionReturn(0); 384d5b485abSSatish Balay } 385d5b485abSSatish Balay /* 386df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 387d5b485abSSatish Balay and return the output 388d5b485abSSatish Balay 389d5b485abSSatish Balay Input: 390d5b485abSSatish Balay C - the matrix 391d5b485abSSatish Balay nrqr - no of messages being processed. 392d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 393d5b485abSSatish Balay 394d5b485abSSatish Balay Output: 395d5b485abSSatish Balay xdata - array of messages to be sent back 396d5b485abSSatish Balay isz1 - size of each message 397d5b485abSSatish Balay 398a8c7a070SBarry Smith For better efficiency perhaps we should malloc separately each xdata[i], 399d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 4000e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of 401d5b485abSSatish Balay memory is used. 402d5b485abSSatish Balay 403d5b485abSSatish Balay */ 404b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 405d5b485abSSatish Balay { 406df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 407d5b485abSSatish Balay Mat A = c->A,B = c->B; 408df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 4096849ba73SBarry Smith PetscErrorCode ierr; 410b24ad042SBarry Smith PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 411b24ad042SBarry Smith PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 412d2a212eaSBarry Smith PetscInt val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 413b24ad042SBarry Smith PetscInt *rbuf_i,kmax,rbuf_0; 414f1af5d2fSBarry Smith PetscBT xtable; 415d5b485abSSatish Balay 4163a40ed3dSBarry Smith PetscFunctionBegin; 417df36731bSBarry Smith Mbs = c->Mbs; 418899cda47SBarry Smith rstart = c->rstartbs; 419899cda47SBarry Smith cstart = c->cstartbs; 420d5b485abSSatish Balay ai = a->i; 421df36731bSBarry Smith aj = a->j; 422d5b485abSSatish Balay bi = b->i; 423df36731bSBarry Smith bj = b->j; 424d5b485abSSatish Balay garray = c->garray; 425d5b485abSSatish Balay 426d5b485abSSatish Balay 427d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 428d5b485abSSatish Balay rbuf_i = rbuf[i]; 429d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 430d5b485abSSatish Balay ct += rbuf_0; 43126fbe8dcSKarl Rupp for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j]; 432d5b485abSSatish Balay } 433d5b485abSSatish Balay 434701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 435701b8814SBarry Smith else max1 = 1; 436d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 437785e854fSJed Brown ierr = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr); 438d5b485abSSatish Balay ++no_malloc; 43953b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr); 440580bdb30SBarry Smith ierr = PetscArrayzero(isz1,nrqr);CHKERRQ(ierr); 441d5b485abSSatish Balay 442d5b485abSSatish Balay ct3 = 0; 443d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 444d5b485abSSatish Balay rbuf_i = rbuf[i]; 445d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 446d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 447d5b485abSSatish Balay ct2 = ct1; 448d5b485abSSatish Balay ct3 += ct1; 449d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 4506831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 451d5b485abSSatish Balay oct2 = ct2; 452d5b485abSSatish Balay kmax = rbuf_i[2*j]; 453d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 454d5b485abSSatish Balay row = rbuf_i[ct1]; 455f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 456d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 457b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 458785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 459580bdb30SBarry Smith ierr = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr); 460606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 461d5b485abSSatish Balay xdata[0] = tmp; 462d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 46326fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 464d5b485abSSatish Balay } 465d5b485abSSatish Balay xdata[i][ct2++] = row; 466d5b485abSSatish Balay ct3++; 467d5b485abSSatish Balay } 468d5b485abSSatish Balay } 469d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 470d5b485abSSatish Balay row = xdata[i][k] - rstart; 471d5b485abSSatish Balay start = ai[row]; 472d5b485abSSatish Balay end = ai[row+1]; 473d5b485abSSatish Balay for (l=start; l<end; l++) { 474df36731bSBarry Smith val = aj[l] + cstart; 475f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 476d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 477b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 478785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 479580bdb30SBarry Smith ierr = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr); 480606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 481d5b485abSSatish Balay xdata[0] = tmp; 482d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 48326fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 484d5b485abSSatish Balay } 485d5b485abSSatish Balay xdata[i][ct2++] = val; 486d5b485abSSatish Balay ct3++; 487d5b485abSSatish Balay } 488d5b485abSSatish Balay } 489d5b485abSSatish Balay start = bi[row]; 490d5b485abSSatish Balay end = bi[row+1]; 491d5b485abSSatish Balay for (l=start; l<end; l++) { 492df36731bSBarry Smith val = garray[bj[l]]; 493f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 494d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 495b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 496785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 497580bdb30SBarry Smith ierr = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr); 498606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 499d5b485abSSatish Balay xdata[0] = tmp; 500d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 50126fbe8dcSKarl Rupp for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 502d5b485abSSatish Balay } 503d5b485abSSatish Balay xdata[i][ct2++] = val; 504d5b485abSSatish Balay ct3++; 505d5b485abSSatish Balay } 506d5b485abSSatish Balay } 507d5b485abSSatish Balay } 508d5b485abSSatish Balay /* Update the header*/ 509d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 510d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 511d5b485abSSatish Balay } 512d5b485abSSatish Balay xdata[i][0] = rbuf_0; 513d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 514d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 515d5b485abSSatish Balay } 51694bacf5dSBarry Smith ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr); 5171e2582c4SBarry Smith ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr); 5183a40ed3dSBarry Smith PetscFunctionReturn(0); 519d5b485abSSatish Balay } 520d5b485abSSatish Balay 5217dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 522d5b485abSSatish Balay { 52317df9f7cSHong Zhang IS *isrow_block,*iscol_block; 524cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 5256849ba73SBarry Smith PetscErrorCode ierr; 526121971b2SHong Zhang PetscInt nmax,nstages_local,nstages,i,pos,max_no,N=C->cmap->N,bs=C->rmap->bs; 5275dc98d6aSHong Zhang Mat_SeqBAIJ *subc; 5285c39f6d9SHong Zhang Mat_SubSppt *smat; 529a2feddc5SSatish Balay 5303a40ed3dSBarry Smith PetscFunctionBegin; 53129dcf524SDmitry Karpeev /* The compression and expansion should be avoided. Doesn't point 53229dcf524SDmitry Karpeev out errors, might change the indices, hence buggey */ 53317df9f7cSHong Zhang ierr = PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block);CHKERRQ(ierr); 53417df9f7cSHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block);CHKERRQ(ierr); 53517df9f7cSHong Zhang ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block);CHKERRQ(ierr); 536cf4f527aSSatish Balay 537cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 5384698041cSHong Zhang if (!C->cmap->N) nmax=20*1000000/sizeof(PetscInt); 5394698041cSHong Zhang else nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 540329f5518SBarry Smith if (!nmax) nmax = 1; 5415dc98d6aSHong Zhang 5425dc98d6aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 5434698041cSHong Zhang nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */ 544cf4f527aSSatish Balay 545653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 546b2566f29SBarry Smith ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 5475dc98d6aSHong Zhang 5484698041cSHong Zhang /* Allocate memory to hold all the submatrices and dummy submatrices */ 5494698041cSHong Zhang ierr = PetscCalloc1(ismax+nstages,submat);CHKERRQ(ierr); 5504698041cSHong Zhang } else { /* MAT_REUSE_MATRIX */ 5514698041cSHong Zhang if (ismax) { 5525dc98d6aSHong Zhang subc = (Mat_SeqBAIJ*)((*submat)[0]->data); 5535dc98d6aSHong Zhang smat = subc->submatis1; 5544698041cSHong Zhang } else { /* (*submat)[0] is a dummy matrix */ 5555c39f6d9SHong Zhang smat = (Mat_SubSppt*)(*submat)[0]->data; 5564698041cSHong Zhang } 557*071fcb05SBarry Smith if (!smat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat"); 5585dc98d6aSHong Zhang nstages = smat->nstages; 5595dc98d6aSHong Zhang } 5605dc98d6aSHong Zhang 561cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 562cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 563cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 564cf4f527aSSatish Balay else max_no = ismax-pos; 565a42ab0d6SHong Zhang 5667dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,*submat+pos);CHKERRQ(ierr); 5674698041cSHong Zhang if (!max_no && scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */ 5685c39f6d9SHong Zhang smat = (Mat_SubSppt*)(*submat)[pos]->data; 5694698041cSHong Zhang smat->nstages = nstages; 5704698041cSHong Zhang } 571cf4f527aSSatish Balay pos += max_no; 572cf4f527aSSatish Balay } 57336f4e84dSSatish Balay 5745dc98d6aSHong Zhang if (scall == MAT_INITIAL_MATRIX && ismax) { 5755dc98d6aSHong Zhang /* save nstages for reuse */ 5765dc98d6aSHong Zhang subc = (Mat_SeqBAIJ*)((*submat)[0]->data); 5775dc98d6aSHong Zhang smat = subc->submatis1; 5785dc98d6aSHong Zhang smat->nstages = nstages; 5795dc98d6aSHong Zhang } 5805dc98d6aSHong Zhang 58136f4e84dSSatish Balay for (i=0; i<ismax; i++) { 58217df9f7cSHong Zhang ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr); 58317df9f7cSHong Zhang ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr); 58436f4e84dSSatish Balay } 58517df9f7cSHong Zhang ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr); 5863a40ed3dSBarry Smith PetscFunctionReturn(0); 587a2feddc5SSatish Balay } 588a2feddc5SSatish Balay 589233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 590e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 591233c2ff6SSatish Balay { 592e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 5934dc2109aSBarry Smith PetscMPIInt fproc; 5944dc2109aSBarry Smith PetscErrorCode ierr; 595233c2ff6SSatish Balay 596233c2ff6SSatish Balay PetscFunctionBegin; 5974dc2109aSBarry Smith ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr); 59823ce1328SBarry Smith if (fproc > size) fproc = size; 599e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 600e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 601233c2ff6SSatish Balay else fproc++; 602233c2ff6SSatish Balay } 603e005ede5SBarry Smith *rank = fproc; 604233c2ff6SSatish Balay PetscFunctionReturn(0); 605233c2ff6SSatish Balay } 606233c2ff6SSatish Balay #endif 607233c2ff6SSatish Balay 608a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 609b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 6107dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats) 61117df9f7cSHong Zhang { 61217df9f7cSHong Zhang Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 61317df9f7cSHong Zhang Mat A = c->A; 61417df9f7cSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc; 61517df9f7cSHong Zhang const PetscInt **icol,**irow; 61617df9f7cSHong Zhang PetscInt *nrow,*ncol,start; 61717df9f7cSHong Zhang PetscErrorCode ierr; 61817df9f7cSHong Zhang PetscMPIInt rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr; 619ddb3d6bcSHong Zhang PetscInt **sbuf1,**sbuf2,*sbuf2_i,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1; 62017df9f7cSHong Zhang PetscInt nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol; 62117df9f7cSHong Zhang PetscInt **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2; 62217df9f7cSHong Zhang PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 62317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 62417df9f7cSHong Zhang PetscTable *cmap,cmap_i=NULL,*rmap,rmap_i; 62517df9f7cSHong Zhang #else 62617df9f7cSHong Zhang PetscInt **cmap,*cmap_i=NULL,**rmap,*rmap_i; 62717df9f7cSHong Zhang #endif 62896cff833SHong Zhang const PetscInt *irow_i,*icol_i; 62917df9f7cSHong Zhang PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 63017df9f7cSHong Zhang MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 63117df9f7cSHong Zhang MPI_Request *r_waits4,*s_waits3,*s_waits4; 63217df9f7cSHong Zhang MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 63317df9f7cSHong Zhang MPI_Status *r_status3,*r_status4,*s_status4; 63417df9f7cSHong Zhang MPI_Comm comm; 63564f5bb2dSSatish Balay PetscScalar **rbuf4,*rbuf4_i=NULL,**sbuf_aa,*vals,*mat_a=NULL,*imat_a=NULL,*sbuf_aa_i; 63617df9f7cSHong Zhang PetscMPIInt *onodes1,*olengths1,end; 63780d31651SHong Zhang PetscInt **row2proc,*row2proc_i,*imat_ilen,*imat_j,*imat_i; 6385c39f6d9SHong Zhang Mat_SubSppt *smat_i; 63917df9f7cSHong Zhang PetscBool *issorted,colflag,iscsorted=PETSC_TRUE; 64017df9f7cSHong Zhang PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen; 641a42ab0d6SHong Zhang PetscInt bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs; 642a42ab0d6SHong Zhang PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 643a42ab0d6SHong Zhang PetscInt nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB; 644afb3cbd8SHong Zhang PetscScalar *vworkA=NULL,*vworkB=NULL,*a_a = a->a,*b_a = b->a; 645a42ab0d6SHong Zhang PetscInt cstart = c->cstartbs,*bmap = c->garray; 64680d31651SHong Zhang PetscBool *allrows,*allcolumns; 647a42ab0d6SHong Zhang 64817df9f7cSHong Zhang PetscFunctionBegin; 64917df9f7cSHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 65017df9f7cSHong Zhang size = c->size; 65117df9f7cSHong Zhang rank = c->rank; 65217df9f7cSHong Zhang 6534698041cSHong Zhang ierr = PetscMalloc5(ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax+1,&allcolumns,ismax,&allrows);CHKERRQ(ierr); 654f489ac74SBarry Smith ierr = PetscMalloc5(ismax,(PetscInt***)&irow,ismax,(PetscInt***)&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr); 65517df9f7cSHong Zhang 65617df9f7cSHong Zhang for (i=0; i<ismax; i++) { 65717df9f7cSHong Zhang ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr); 65880d31651SHong Zhang if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */ 65917df9f7cSHong Zhang ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr); 66017df9f7cSHong Zhang 661ec3c739cSHong Zhang /* Check for special case: allcolumns */ 66217df9f7cSHong Zhang ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 66317df9f7cSHong Zhang ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 6645dd0c08cSHong Zhang 6651abc2f7bSHong Zhang if (colflag && ncol[i] == c->Nbs) { 6661abc2f7bSHong Zhang allcolumns[i] = PETSC_TRUE; 6671abc2f7bSHong Zhang icol[i] = NULL; 6681abc2f7bSHong Zhang } else { 66917df9f7cSHong Zhang allcolumns[i] = PETSC_FALSE; 67017df9f7cSHong Zhang ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 6711abc2f7bSHong Zhang } 672ec3c739cSHong Zhang 673ec3c739cSHong Zhang /* Check for special case: allrows */ 674ec3c739cSHong Zhang ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 675ec3c739cSHong Zhang ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 676ec3c739cSHong Zhang if (colflag && nrow[i] == c->Mbs) { 677ec3c739cSHong Zhang allrows[i] = PETSC_TRUE; 678ec3c739cSHong Zhang irow[i] = NULL; 679ec3c739cSHong Zhang } else { 680ec3c739cSHong Zhang allrows[i] = PETSC_FALSE; 681ec3c739cSHong Zhang ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 682ec3c739cSHong Zhang } 68317df9f7cSHong Zhang } 68417df9f7cSHong Zhang 68517df9f7cSHong Zhang if (scall == MAT_REUSE_MATRIX) { 68617df9f7cSHong Zhang /* Assumes new rows are same length as the old rows */ 68717df9f7cSHong Zhang for (i=0; i<ismax; i++) { 68817df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)(submats[i]->data); 6899e686edcSHong Zhang if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 69017df9f7cSHong Zhang 69117df9f7cSHong Zhang /* Initial matrix as if empty */ 692580bdb30SBarry Smith ierr = PetscArrayzero(subc->ilen,subc->mbs);CHKERRQ(ierr); 69317df9f7cSHong Zhang 69417df9f7cSHong Zhang /* Initial matrix as if empty */ 69517df9f7cSHong Zhang submats[i]->factortype = C->factortype; 69617df9f7cSHong Zhang 69717df9f7cSHong Zhang smat_i = subc->submatis1; 69817df9f7cSHong Zhang 69917df9f7cSHong Zhang nrqs = smat_i->nrqs; 70017df9f7cSHong Zhang nrqr = smat_i->nrqr; 70117df9f7cSHong Zhang rbuf1 = smat_i->rbuf1; 70217df9f7cSHong Zhang rbuf2 = smat_i->rbuf2; 70317df9f7cSHong Zhang rbuf3 = smat_i->rbuf3; 70417df9f7cSHong Zhang req_source2 = smat_i->req_source2; 70517df9f7cSHong Zhang 70617df9f7cSHong Zhang sbuf1 = smat_i->sbuf1; 70717df9f7cSHong Zhang sbuf2 = smat_i->sbuf2; 70817df9f7cSHong Zhang ptr = smat_i->ptr; 70917df9f7cSHong Zhang tmp = smat_i->tmp; 71017df9f7cSHong Zhang ctr = smat_i->ctr; 71117df9f7cSHong Zhang 71217df9f7cSHong Zhang pa = smat_i->pa; 71317df9f7cSHong Zhang req_size = smat_i->req_size; 71417df9f7cSHong Zhang req_source1 = smat_i->req_source1; 71517df9f7cSHong Zhang 71617df9f7cSHong Zhang allcolumns[i] = smat_i->allcolumns; 717ec3c739cSHong Zhang allrows[i] = smat_i->allrows; 71817df9f7cSHong Zhang row2proc[i] = smat_i->row2proc; 71917df9f7cSHong Zhang rmap[i] = smat_i->rmap; 72017df9f7cSHong Zhang cmap[i] = smat_i->cmap; 72117df9f7cSHong Zhang } 7224698041cSHong Zhang 7234698041cSHong Zhang if (!ismax){ /* Get dummy submatrices and retrieve struct submatis1 */ 7244698041cSHong Zhang if (!submats[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"submats are null, cannot reuse"); 7255c39f6d9SHong Zhang smat_i = (Mat_SubSppt*)submats[0]->data; 7264698041cSHong Zhang 7274698041cSHong Zhang nrqs = smat_i->nrqs; 7284698041cSHong Zhang nrqr = smat_i->nrqr; 7294698041cSHong Zhang rbuf1 = smat_i->rbuf1; 7304698041cSHong Zhang rbuf2 = smat_i->rbuf2; 7314698041cSHong Zhang rbuf3 = smat_i->rbuf3; 7324698041cSHong Zhang req_source2 = smat_i->req_source2; 7334698041cSHong Zhang 7344698041cSHong Zhang sbuf1 = smat_i->sbuf1; 7354698041cSHong Zhang sbuf2 = smat_i->sbuf2; 7364698041cSHong Zhang ptr = smat_i->ptr; 7374698041cSHong Zhang tmp = smat_i->tmp; 7384698041cSHong Zhang ctr = smat_i->ctr; 7394698041cSHong Zhang 7404698041cSHong Zhang pa = smat_i->pa; 7414698041cSHong Zhang req_size = smat_i->req_size; 7424698041cSHong Zhang req_source1 = smat_i->req_source1; 7434698041cSHong Zhang 744c67fe082SHong Zhang allcolumns[0] = PETSC_FALSE; 7454698041cSHong Zhang } 74617df9f7cSHong Zhang } else { /* scall == MAT_INITIAL_MATRIX */ 74717df9f7cSHong Zhang /* Get some new tags to keep the communication clean */ 74817df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 74917df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 75017df9f7cSHong Zhang 75117df9f7cSHong Zhang /* evaluate communication - mesg to who, length of mesg, and buffer space 75217df9f7cSHong Zhang required. Based on this, buffers are allocated, and data copied into them*/ 75317df9f7cSHong Zhang ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); /* mesg size, initialize work vectors */ 75417df9f7cSHong Zhang 75517df9f7cSHong Zhang for (i=0; i<ismax; i++) { 75617df9f7cSHong Zhang jmax = nrow[i]; 75717df9f7cSHong Zhang irow_i = irow[i]; 75817df9f7cSHong Zhang 75917df9f7cSHong Zhang ierr = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr); 76017df9f7cSHong Zhang row2proc[i] = row2proc_i; 76117df9f7cSHong Zhang 76217df9f7cSHong Zhang if (issorted[i]) proc = 0; 76317df9f7cSHong Zhang for (j=0; j<jmax; j++) { 76417df9f7cSHong Zhang if (!issorted[i]) proc = 0; 765ec3c739cSHong Zhang if (allrows[i]) row = j; 766ec3c739cSHong Zhang else row = irow_i[j]; 767ec3c739cSHong Zhang 768a42ab0d6SHong Zhang while (row >= c->rangebs[proc+1]) proc++; 76917df9f7cSHong Zhang w4[proc]++; 77017df9f7cSHong Zhang row2proc_i[j] = proc; /* map row index to proc */ 77117df9f7cSHong Zhang } 77217df9f7cSHong Zhang for (j=0; j<size; j++) { 77317df9f7cSHong Zhang if (w4[j]) { w1[j] += w4[j]; w3[j]++; w4[j] = 0;} 77417df9f7cSHong Zhang } 77517df9f7cSHong Zhang } 77617df9f7cSHong Zhang 77717df9f7cSHong Zhang nrqs = 0; /* no of outgoing messages */ 77817df9f7cSHong Zhang msz = 0; /* total mesg length (for all procs) */ 77917df9f7cSHong Zhang w1[rank] = 0; /* no mesg sent to self */ 78017df9f7cSHong Zhang w3[rank] = 0; 78117df9f7cSHong Zhang for (i=0; i<size; i++) { 78217df9f7cSHong Zhang if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 78317df9f7cSHong Zhang } 78417df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 78517df9f7cSHong Zhang for (i=0,j=0; i<size; i++) { 78617df9f7cSHong Zhang if (w1[i]) { pa[j] = i; j++; } 78717df9f7cSHong Zhang } 78817df9f7cSHong Zhang 78917df9f7cSHong Zhang /* Each message would have a header = 1 + 2*(no of IS) + data */ 79017df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 79117df9f7cSHong Zhang j = pa[i]; 79217df9f7cSHong Zhang w1[j] += w2[j] + 2* w3[j]; 79317df9f7cSHong Zhang msz += w1[j]; 79417df9f7cSHong Zhang } 79517df9f7cSHong Zhang ierr = PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);CHKERRQ(ierr); 79617df9f7cSHong Zhang 79717df9f7cSHong Zhang /* Determine the number of messages to expect, their lengths, from from-ids */ 79817df9f7cSHong Zhang ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 79917df9f7cSHong Zhang ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 80017df9f7cSHong Zhang 80117df9f7cSHong Zhang /* Now post the Irecvs corresponding to these messages */ 80217df9f7cSHong Zhang tag0 = ((PetscObject)C)->tag; 80317df9f7cSHong Zhang ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 80417df9f7cSHong Zhang 80517df9f7cSHong Zhang ierr = PetscFree(onodes1);CHKERRQ(ierr); 80617df9f7cSHong Zhang ierr = PetscFree(olengths1);CHKERRQ(ierr); 80717df9f7cSHong Zhang 80817df9f7cSHong Zhang /* Allocate Memory for outgoing messages */ 80917df9f7cSHong Zhang ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 810580bdb30SBarry Smith ierr = PetscArrayzero(sbuf1,size);CHKERRQ(ierr); 811580bdb30SBarry Smith ierr = PetscArrayzero(ptr,size);CHKERRQ(ierr); 81217df9f7cSHong Zhang 81317df9f7cSHong Zhang { 81417df9f7cSHong Zhang PetscInt *iptr = tmp; 81517df9f7cSHong Zhang k = 0; 81617df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 81717df9f7cSHong Zhang j = pa[i]; 81817df9f7cSHong Zhang iptr += k; 81917df9f7cSHong Zhang sbuf1[j] = iptr; 82017df9f7cSHong Zhang k = w1[j]; 82117df9f7cSHong Zhang } 82217df9f7cSHong Zhang } 82317df9f7cSHong Zhang 82417df9f7cSHong Zhang /* Form the outgoing messages. Initialize the header space */ 82517df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 82617df9f7cSHong Zhang j = pa[i]; 82717df9f7cSHong Zhang sbuf1[j][0] = 0; 828580bdb30SBarry Smith ierr = PetscArrayzero(sbuf1[j]+1,2*w3[j]);CHKERRQ(ierr); 82917df9f7cSHong Zhang ptr[j] = sbuf1[j] + 2*w3[j] + 1; 83017df9f7cSHong Zhang } 83117df9f7cSHong Zhang 83217df9f7cSHong Zhang /* Parse the isrow and copy data into outbuf */ 83317df9f7cSHong Zhang for (i=0; i<ismax; i++) { 83417df9f7cSHong Zhang row2proc_i = row2proc[i]; 835580bdb30SBarry Smith ierr = PetscArrayzero(ctr,size);CHKERRQ(ierr); 83617df9f7cSHong Zhang irow_i = irow[i]; 83717df9f7cSHong Zhang jmax = nrow[i]; 83817df9f7cSHong Zhang for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 83917df9f7cSHong Zhang proc = row2proc_i[j]; 840ec3c739cSHong Zhang if (allrows[i]) row = j; 841ec3c739cSHong Zhang else row = irow_i[j]; 842ec3c739cSHong Zhang 84317df9f7cSHong Zhang if (proc != rank) { /* copy to the outgoing buf*/ 84417df9f7cSHong Zhang ctr[proc]++; 845ec3c739cSHong Zhang *ptr[proc] = row; 84617df9f7cSHong Zhang ptr[proc]++; 84717df9f7cSHong Zhang } 84817df9f7cSHong Zhang } 84917df9f7cSHong Zhang /* Update the headers for the current IS */ 85017df9f7cSHong Zhang for (j=0; j<size; j++) { /* Can Optimise this loop too */ 85117df9f7cSHong Zhang if ((ctr_j = ctr[j])) { 85217df9f7cSHong Zhang sbuf1_j = sbuf1[j]; 85317df9f7cSHong Zhang k = ++sbuf1_j[0]; 85417df9f7cSHong Zhang sbuf1_j[2*k] = ctr_j; 85517df9f7cSHong Zhang sbuf1_j[2*k-1] = i; 85617df9f7cSHong Zhang } 85717df9f7cSHong Zhang } 85817df9f7cSHong Zhang } 85917df9f7cSHong Zhang 86017df9f7cSHong Zhang /* Now post the sends */ 86117df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 86217df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 86317df9f7cSHong Zhang j = pa[i]; 86417df9f7cSHong Zhang ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 86517df9f7cSHong Zhang } 86617df9f7cSHong Zhang 86717df9f7cSHong Zhang /* Post Receives to capture the buffer size */ 86817df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 86917df9f7cSHong Zhang ierr = PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);CHKERRQ(ierr); 87017df9f7cSHong Zhang rbuf2[0] = tmp + msz; 87117df9f7cSHong Zhang for (i=1; i<nrqs; ++i) { 87217df9f7cSHong Zhang rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 87317df9f7cSHong Zhang } 87417df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 87517df9f7cSHong Zhang j = pa[i]; 87617df9f7cSHong Zhang ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRQ(ierr); 87717df9f7cSHong Zhang } 87817df9f7cSHong Zhang 87917df9f7cSHong Zhang /* Send to other procs the buf size they should allocate */ 88017df9f7cSHong Zhang /* Receive messages*/ 88117df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 88217df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 88317df9f7cSHong Zhang ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr); 88417df9f7cSHong Zhang 88517df9f7cSHong Zhang ierr = MPI_Waitall(nrqr,r_waits1,r_status1);CHKERRQ(ierr); 88617df9f7cSHong Zhang for (i=0; i<nrqr; ++i) { 88717df9f7cSHong Zhang req_size[i] = 0; 88817df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 88917df9f7cSHong Zhang start = 2*rbuf1_i[0] + 1; 89017df9f7cSHong Zhang ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 89117df9f7cSHong Zhang ierr = PetscMalloc1(end+1,&sbuf2[i]);CHKERRQ(ierr); 89217df9f7cSHong Zhang sbuf2_i = sbuf2[i]; 89317df9f7cSHong Zhang for (j=start; j<end; j++) { 894ddb3d6bcSHong Zhang row = rbuf1_i[j] - rstart; 895ddb3d6bcSHong Zhang ncols = a_i[row+1] - a_i[row] + b_i[row+1] - b_i[row]; 89617df9f7cSHong Zhang sbuf2_i[j] = ncols; 89717df9f7cSHong Zhang req_size[i] += ncols; 89817df9f7cSHong Zhang } 89917df9f7cSHong Zhang req_source1[i] = r_status1[i].MPI_SOURCE; 90017df9f7cSHong Zhang /* form the header */ 90117df9f7cSHong Zhang sbuf2_i[0] = req_size[i]; 90217df9f7cSHong Zhang for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 90317df9f7cSHong Zhang 90417df9f7cSHong Zhang ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRQ(ierr); 90517df9f7cSHong Zhang } 906ddb3d6bcSHong Zhang 90717df9f7cSHong Zhang ierr = PetscFree(r_status1);CHKERRQ(ierr); 90817df9f7cSHong Zhang ierr = PetscFree(r_waits1);CHKERRQ(ierr); 90917df9f7cSHong Zhang ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 91017df9f7cSHong Zhang 91117df9f7cSHong Zhang /* Receive messages*/ 91217df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 91317df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 91417df9f7cSHong Zhang 91517df9f7cSHong Zhang ierr = MPI_Waitall(nrqs,r_waits2,r_status2);CHKERRQ(ierr); 91617df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 91717df9f7cSHong Zhang ierr = PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);CHKERRQ(ierr); 91817df9f7cSHong Zhang req_source2[i] = r_status2[i].MPI_SOURCE; 91917df9f7cSHong Zhang ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRQ(ierr); 92017df9f7cSHong Zhang } 92117df9f7cSHong Zhang ierr = PetscFree(r_status2);CHKERRQ(ierr); 92217df9f7cSHong Zhang ierr = PetscFree(r_waits2);CHKERRQ(ierr); 92317df9f7cSHong Zhang 92417df9f7cSHong Zhang /* Wait on sends1 and sends2 */ 92517df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 92617df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 92717df9f7cSHong Zhang 92817df9f7cSHong Zhang if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 92917df9f7cSHong Zhang if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 93017df9f7cSHong Zhang ierr = PetscFree(s_status1);CHKERRQ(ierr); 93117df9f7cSHong Zhang ierr = PetscFree(s_status2);CHKERRQ(ierr); 93217df9f7cSHong Zhang ierr = PetscFree(s_waits1);CHKERRQ(ierr); 93317df9f7cSHong Zhang ierr = PetscFree(s_waits2);CHKERRQ(ierr); 93417df9f7cSHong Zhang 93517df9f7cSHong Zhang /* Now allocate sending buffers for a->j, and send them off */ 93617df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 93717df9f7cSHong Zhang for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 93817df9f7cSHong Zhang ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 93917df9f7cSHong Zhang for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 94017df9f7cSHong Zhang 94117df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 94217df9f7cSHong Zhang { 94317df9f7cSHong Zhang 94417df9f7cSHong Zhang for (i=0; i<nrqr; i++) { 94517df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 94617df9f7cSHong Zhang sbuf_aj_i = sbuf_aj[i]; 94717df9f7cSHong Zhang ct1 = 2*rbuf1_i[0] + 1; 94817df9f7cSHong Zhang ct2 = 0; 94917df9f7cSHong Zhang for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 95017df9f7cSHong Zhang kmax = rbuf1[i][2*j]; 95117df9f7cSHong Zhang for (k=0; k<kmax; k++,ct1++) { 95217df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 95317df9f7cSHong Zhang nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 95417df9f7cSHong Zhang ncols = nzA + nzB; 95517df9f7cSHong Zhang cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 95617df9f7cSHong Zhang 95717df9f7cSHong Zhang /* load the column indices for this row into cols */ 95817df9f7cSHong Zhang cols = sbuf_aj_i + ct2; 95917df9f7cSHong Zhang for (l=0; l<nzB; l++) { 960a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 961a42ab0d6SHong Zhang else break; 96217df9f7cSHong Zhang } 963a42ab0d6SHong Zhang imark = l; 964a42ab0d6SHong Zhang for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];} 965a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 96617df9f7cSHong Zhang ct2 += ncols; 96717df9f7cSHong Zhang } 96817df9f7cSHong Zhang } 96917df9f7cSHong Zhang ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRQ(ierr); 97017df9f7cSHong Zhang } 97117df9f7cSHong Zhang } 97217df9f7cSHong Zhang ierr = PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);CHKERRQ(ierr); 97317df9f7cSHong Zhang 97417df9f7cSHong Zhang /* create col map: global col of C -> local col of submatrices */ 97517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 97617df9f7cSHong Zhang for (i=0; i<ismax; i++) { 97717df9f7cSHong Zhang if (!allcolumns[i]) { 978a42ab0d6SHong Zhang ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 97917df9f7cSHong Zhang 98017df9f7cSHong Zhang jmax = ncol[i]; 98117df9f7cSHong Zhang icol_i = icol[i]; 98217df9f7cSHong Zhang cmap_i = cmap[i]; 98317df9f7cSHong Zhang for (j=0; j<jmax; j++) { 98417df9f7cSHong Zhang ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 98517df9f7cSHong Zhang } 98617df9f7cSHong Zhang } else cmap[i] = NULL; 98717df9f7cSHong Zhang } 98817df9f7cSHong Zhang #else 98917df9f7cSHong Zhang for (i=0; i<ismax; i++) { 99017df9f7cSHong Zhang if (!allcolumns[i]) { 991a42ab0d6SHong Zhang ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 99217df9f7cSHong Zhang jmax = ncol[i]; 99317df9f7cSHong Zhang icol_i = icol[i]; 99417df9f7cSHong Zhang cmap_i = cmap[i]; 995a42ab0d6SHong Zhang for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1; 99617df9f7cSHong Zhang } else cmap[i] = NULL; 99717df9f7cSHong Zhang } 99817df9f7cSHong Zhang #endif 99917df9f7cSHong Zhang 100017df9f7cSHong Zhang /* Create lens which is required for MatCreate... */ 100117df9f7cSHong Zhang for (i=0,j=0; i<ismax; i++) j += nrow[i]; 100217df9f7cSHong Zhang ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr); 100317df9f7cSHong Zhang 100417df9f7cSHong Zhang if (ismax) { 100517df9f7cSHong Zhang ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr); 100617df9f7cSHong Zhang } 100717df9f7cSHong Zhang for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 100817df9f7cSHong Zhang 100917df9f7cSHong Zhang /* Update lens from local data */ 101017df9f7cSHong Zhang for (i=0; i<ismax; i++) { 101117df9f7cSHong Zhang row2proc_i = row2proc[i]; 101217df9f7cSHong Zhang jmax = nrow[i]; 101317df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 101417df9f7cSHong Zhang irow_i = irow[i]; 101517df9f7cSHong Zhang lens_i = lens[i]; 101617df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1017ec3c739cSHong Zhang if (allrows[i]) row = j; 1018ec3c739cSHong Zhang else row = irow_i[j]; /* global blocked row of C */ 1019ec3c739cSHong Zhang 102017df9f7cSHong Zhang proc = row2proc_i[j]; 102117df9f7cSHong Zhang if (proc == rank) { 1022a42ab0d6SHong Zhang /* Get indices from matA and then from matB */ 102317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1024a42ab0d6SHong Zhang PetscInt tt; 1025a42ab0d6SHong Zhang #endif 1026a42ab0d6SHong Zhang row = row - rstart; 1027a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1028a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 1029a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 1030a42ab0d6SHong Zhang cworkB = b_j + b_i[row]; 1031a42ab0d6SHong Zhang 1032a42ab0d6SHong Zhang if (!allcolumns[i]) { 1033a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1034a42ab0d6SHong Zhang for (k=0; k<nzA; k++) { 1035a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1036a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1037a42ab0d6SHong Zhang } 1038a42ab0d6SHong Zhang for (k=0; k<nzB; k++) { 1039a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1040a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1041a42ab0d6SHong Zhang } 1042a42ab0d6SHong Zhang 1043a42ab0d6SHong Zhang #else 1044a42ab0d6SHong Zhang for (k=0; k<nzA; k++) { 1045a42ab0d6SHong Zhang if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1046a42ab0d6SHong Zhang } 1047a42ab0d6SHong Zhang for (k=0; k<nzB; k++) { 1048a42ab0d6SHong Zhang if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1049a42ab0d6SHong Zhang } 1050a42ab0d6SHong Zhang #endif 1051a42ab0d6SHong Zhang } else { /* allcolumns */ 1052a42ab0d6SHong Zhang lens_i[j] = nzA + nzB; 1053a42ab0d6SHong Zhang } 105417df9f7cSHong Zhang } 105517df9f7cSHong Zhang } 105617df9f7cSHong Zhang } 105717df9f7cSHong Zhang 105817df9f7cSHong Zhang /* Create row map: global row of C -> local row of submatrices */ 105917df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1060059f6b74SHong Zhang if (!allrows[i]) { 1061059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE) 1062a42ab0d6SHong Zhang ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 106317df9f7cSHong Zhang irow_i = irow[i]; 106417df9f7cSHong Zhang jmax = nrow[i]; 106517df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1066ec3c739cSHong Zhang if (allrows[i]) { 1067ec3c739cSHong Zhang ierr = PetscTableAdd(rmap[i],j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1068ec3c739cSHong Zhang } else { 106917df9f7cSHong Zhang ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 107017df9f7cSHong Zhang } 107117df9f7cSHong Zhang } 107217df9f7cSHong Zhang #else 1073a42ab0d6SHong Zhang ierr = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr); 107417df9f7cSHong Zhang rmap_i = rmap[i]; 107517df9f7cSHong Zhang irow_i = irow[i]; 107617df9f7cSHong Zhang jmax = nrow[i]; 107717df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1078ec3c739cSHong Zhang if (allrows[i]) rmap_i[j] = j; 1079ec3c739cSHong Zhang else rmap_i[irow_i[j]] = j; 108017df9f7cSHong Zhang } 108117df9f7cSHong Zhang #endif 1082059f6b74SHong Zhang } else rmap[i] = NULL; 1083059f6b74SHong Zhang } 108417df9f7cSHong Zhang 108517df9f7cSHong Zhang /* Update lens from offproc data */ 108617df9f7cSHong Zhang { 108717df9f7cSHong Zhang PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 108817df9f7cSHong Zhang 108917df9f7cSHong Zhang ierr = MPI_Waitall(nrqs,r_waits3,r_status3);CHKERRQ(ierr); 109017df9f7cSHong Zhang for (tmp2=0; tmp2<nrqs; tmp2++) { 109117df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 109217df9f7cSHong Zhang jmax = sbuf1_i[0]; 109317df9f7cSHong Zhang ct1 = 2*jmax+1; 109417df9f7cSHong Zhang ct2 = 0; 109517df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 109617df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 109717df9f7cSHong Zhang for (j=1; j<=jmax; j++) { 109817df9f7cSHong Zhang is_no = sbuf1_i[2*j-1]; 109917df9f7cSHong Zhang max1 = sbuf1_i[2*j]; 110017df9f7cSHong Zhang lens_i = lens[is_no]; 110117df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 110217df9f7cSHong Zhang rmap_i = rmap[is_no]; 110317df9f7cSHong Zhang for (k=0; k<max1; k++,ct1++) { 1104059f6b74SHong Zhang if (allrows[is_no]) { 1105059f6b74SHong Zhang row = sbuf1_i[ct1]; 1106059f6b74SHong Zhang } else { 110717df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 110817df9f7cSHong Zhang ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 110917df9f7cSHong Zhang row--; 111017df9f7cSHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 111117df9f7cSHong Zhang #else 111217df9f7cSHong Zhang row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 111317df9f7cSHong Zhang #endif 1114059f6b74SHong Zhang } 111517df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 111617df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 111717df9f7cSHong Zhang if (!allcolumns[is_no]) { 111817df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 111917df9f7cSHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 112017df9f7cSHong Zhang #else 112117df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 112217df9f7cSHong Zhang #endif 112317df9f7cSHong Zhang if (tcol) lens_i[row]++; 112417df9f7cSHong Zhang } else { /* allcolumns */ 112517df9f7cSHong Zhang lens_i[row]++; /* lens_i[row] += max2 ? */ 112617df9f7cSHong Zhang } 112717df9f7cSHong Zhang } 112817df9f7cSHong Zhang } 112917df9f7cSHong Zhang } 113017df9f7cSHong Zhang } 113117df9f7cSHong Zhang } 113217df9f7cSHong Zhang ierr = PetscFree(r_waits3);CHKERRQ(ierr); 113317df9f7cSHong Zhang if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 113417df9f7cSHong Zhang ierr = PetscFree2(r_status3,s_status3);CHKERRQ(ierr); 113517df9f7cSHong Zhang ierr = PetscFree(s_waits3);CHKERRQ(ierr); 113617df9f7cSHong Zhang 113717df9f7cSHong Zhang /* Create the submatrices */ 113817df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1139a42ab0d6SHong Zhang PetscInt bs_tmp; 1140a42ab0d6SHong Zhang if (ijonly) bs_tmp = 1; 1141a42ab0d6SHong Zhang else bs_tmp = bs; 114217df9f7cSHong Zhang 114317df9f7cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 1144a42ab0d6SHong Zhang ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 114517df9f7cSHong Zhang 114617df9f7cSHong Zhang ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 1147a42ab0d6SHong Zhang ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 1148a42ab0d6SHong Zhang ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 114917df9f7cSHong Zhang 11505c39f6d9SHong Zhang /* create struct Mat_SubSppt and attached it to submat */ 115117df9f7cSHong Zhang ierr = PetscNew(&smat_i);CHKERRQ(ierr); 115217df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 115317df9f7cSHong Zhang subc->submatis1 = smat_i; 115417df9f7cSHong Zhang 115517df9f7cSHong Zhang smat_i->destroy = submats[i]->ops->destroy; 1156f68bb481SHong Zhang submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ; 115717df9f7cSHong Zhang submats[i]->factortype = C->factortype; 115817df9f7cSHong Zhang 115917df9f7cSHong Zhang smat_i->id = i; 116017df9f7cSHong Zhang smat_i->nrqs = nrqs; 116117df9f7cSHong Zhang smat_i->nrqr = nrqr; 116217df9f7cSHong Zhang smat_i->rbuf1 = rbuf1; 116317df9f7cSHong Zhang smat_i->rbuf2 = rbuf2; 116417df9f7cSHong Zhang smat_i->rbuf3 = rbuf3; 116517df9f7cSHong Zhang smat_i->sbuf2 = sbuf2; 116617df9f7cSHong Zhang smat_i->req_source2 = req_source2; 116717df9f7cSHong Zhang 116817df9f7cSHong Zhang smat_i->sbuf1 = sbuf1; 116917df9f7cSHong Zhang smat_i->ptr = ptr; 117017df9f7cSHong Zhang smat_i->tmp = tmp; 117117df9f7cSHong Zhang smat_i->ctr = ctr; 117217df9f7cSHong Zhang 117317df9f7cSHong Zhang smat_i->pa = pa; 117417df9f7cSHong Zhang smat_i->req_size = req_size; 117517df9f7cSHong Zhang smat_i->req_source1 = req_source1; 117617df9f7cSHong Zhang 117717df9f7cSHong Zhang smat_i->allcolumns = allcolumns[i]; 1178ec3c739cSHong Zhang smat_i->allrows = allrows[i]; 117917df9f7cSHong Zhang smat_i->singleis = PETSC_FALSE; 118017df9f7cSHong Zhang smat_i->row2proc = row2proc[i]; 118117df9f7cSHong Zhang smat_i->rmap = rmap[i]; 118217df9f7cSHong Zhang smat_i->cmap = cmap[i]; 118317df9f7cSHong Zhang } 118417df9f7cSHong Zhang 11854698041cSHong Zhang if (!ismax) { /* Create dummy submats[0] for reuse struct subc */ 11864698041cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&submats[0]);CHKERRQ(ierr); 11874698041cSHong Zhang ierr = MatSetSizes(submats[0],0,0,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 11884698041cSHong Zhang ierr = MatSetType(submats[0],MATDUMMY);CHKERRQ(ierr); 11894698041cSHong Zhang 11905c39f6d9SHong Zhang /* create struct Mat_SubSppt and attached it to submat */ 11914698041cSHong Zhang ierr = PetscNewLog(submats[0],&smat_i);CHKERRQ(ierr); 11924698041cSHong Zhang submats[0]->data = (void*)smat_i; 11934698041cSHong Zhang 11944698041cSHong Zhang smat_i->destroy = submats[0]->ops->destroy; 1195f68bb481SHong Zhang submats[0]->ops->destroy = MatDestroySubMatrix_Dummy; 11964698041cSHong Zhang submats[0]->factortype = C->factortype; 11974698041cSHong Zhang 1198006cef31SHong Zhang smat_i->id = 0; 11994698041cSHong Zhang smat_i->nrqs = nrqs; 12004698041cSHong Zhang smat_i->nrqr = nrqr; 12014698041cSHong Zhang smat_i->rbuf1 = rbuf1; 12024698041cSHong Zhang smat_i->rbuf2 = rbuf2; 12034698041cSHong Zhang smat_i->rbuf3 = rbuf3; 12044698041cSHong Zhang smat_i->sbuf2 = sbuf2; 12054698041cSHong Zhang smat_i->req_source2 = req_source2; 12064698041cSHong Zhang 12074698041cSHong Zhang smat_i->sbuf1 = sbuf1; 12084698041cSHong Zhang smat_i->ptr = ptr; 12094698041cSHong Zhang smat_i->tmp = tmp; 12104698041cSHong Zhang smat_i->ctr = ctr; 12114698041cSHong Zhang 12124698041cSHong Zhang smat_i->pa = pa; 12134698041cSHong Zhang smat_i->req_size = req_size; 12144698041cSHong Zhang smat_i->req_source1 = req_source1; 12154698041cSHong Zhang 1216c67fe082SHong Zhang smat_i->allcolumns = PETSC_FALSE; 12174698041cSHong Zhang smat_i->singleis = PETSC_FALSE; 12184698041cSHong Zhang smat_i->row2proc = NULL; 12194698041cSHong Zhang smat_i->rmap = NULL; 12204698041cSHong Zhang smat_i->cmap = NULL; 12214698041cSHong Zhang } 12224698041cSHong Zhang 12234698041cSHong Zhang 122417df9f7cSHong Zhang if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);} 122517df9f7cSHong Zhang ierr = PetscFree(lens);CHKERRQ(ierr); 122617df9f7cSHong Zhang ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 122717df9f7cSHong Zhang ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 122817df9f7cSHong Zhang 122917df9f7cSHong Zhang } /* endof scall == MAT_INITIAL_MATRIX */ 123017df9f7cSHong Zhang 123117df9f7cSHong Zhang /* Post recv matrix values */ 1232bbc89d27SHong Zhang if (!ijonly) { 123317df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr); 123417df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 123517df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 123617df9f7cSHong Zhang ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 123717df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 123817df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 1239a42ab0d6SHong Zhang ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr); 1240a42ab0d6SHong Zhang ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRQ(ierr); 124117df9f7cSHong Zhang } 124217df9f7cSHong Zhang 124317df9f7cSHong Zhang /* Allocate sending buffers for a->a, and send them off */ 124417df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 124517df9f7cSHong Zhang for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 12469e686edcSHong Zhang 1247a42ab0d6SHong Zhang ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 1248a42ab0d6SHong Zhang for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 124917df9f7cSHong Zhang 125017df9f7cSHong Zhang ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 125117df9f7cSHong Zhang 125217df9f7cSHong Zhang for (i=0; i<nrqr; i++) { 125317df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 125417df9f7cSHong Zhang sbuf_aa_i = sbuf_aa[i]; 125517df9f7cSHong Zhang ct1 = 2*rbuf1_i[0]+1; 125617df9f7cSHong Zhang ct2 = 0; 125717df9f7cSHong Zhang for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 125817df9f7cSHong Zhang kmax = rbuf1_i[2*j]; 125917df9f7cSHong Zhang for (k=0; k<kmax; k++,ct1++) { 126017df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 1261a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1262a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 126317df9f7cSHong Zhang ncols = nzA + nzB; 126417df9f7cSHong Zhang cworkB = b_j + b_i[row]; 1265a42ab0d6SHong Zhang vworkA = a_a + a_i[row]*bs2; 1266a42ab0d6SHong Zhang vworkB = b_a + b_i[row]*bs2; 126717df9f7cSHong Zhang 126817df9f7cSHong Zhang /* load the column values for this row into vals*/ 1269a42ab0d6SHong Zhang vals = sbuf_aa_i+ct2*bs2; 1270a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1271a42ab0d6SHong Zhang if ((bmap[cworkB[l]]) < cstart) { 1272580bdb30SBarry Smith ierr = PetscArraycpy(vals+l*bs2,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1273a42ab0d6SHong Zhang } else break; 1274a42ab0d6SHong Zhang } 1275a42ab0d6SHong Zhang imark = l; 1276a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1277580bdb30SBarry Smith ierr = PetscArraycpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2);CHKERRQ(ierr); 1278a42ab0d6SHong Zhang } 1279a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1280580bdb30SBarry Smith ierr = PetscArraycpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1281a42ab0d6SHong Zhang } 128217df9f7cSHong Zhang 128317df9f7cSHong Zhang ct2 += ncols; 128417df9f7cSHong Zhang } 128517df9f7cSHong Zhang } 1286a42ab0d6SHong Zhang ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRQ(ierr); 128717df9f7cSHong Zhang } 1288bbc89d27SHong Zhang } 128917df9f7cSHong Zhang 129017df9f7cSHong Zhang /* Assemble the matrices */ 129117df9f7cSHong Zhang /* First assemble the local rows */ 129217df9f7cSHong Zhang for (i=0; i<ismax; i++) { 129317df9f7cSHong Zhang row2proc_i = row2proc[i]; 129417df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 129517df9f7cSHong Zhang imat_ilen = subc->ilen; 129617df9f7cSHong Zhang imat_j = subc->j; 129717df9f7cSHong Zhang imat_i = subc->i; 129817df9f7cSHong Zhang imat_a = subc->a; 129917df9f7cSHong Zhang 130017df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 130117df9f7cSHong Zhang rmap_i = rmap[i]; 130217df9f7cSHong Zhang irow_i = irow[i]; 130317df9f7cSHong Zhang jmax = nrow[i]; 130417df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1305ec3c739cSHong Zhang if (allrows[i]) row = j; 1306ec3c739cSHong Zhang else row = irow_i[j]; 130717df9f7cSHong Zhang proc = row2proc_i[j]; 1308a42ab0d6SHong Zhang 130917df9f7cSHong Zhang if (proc == rank) { 1310a42ab0d6SHong Zhang 1311a42ab0d6SHong Zhang row = row - rstart; 1312a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1313a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 1314a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 1315a42ab0d6SHong Zhang cworkB = b_j + b_i[row]; 1316a42ab0d6SHong Zhang if (!ijonly) { 1317a42ab0d6SHong Zhang vworkA = a_a + a_i[row]*bs2; 1318a42ab0d6SHong Zhang vworkB = b_a + b_i[row]*bs2; 1319a42ab0d6SHong Zhang } 1320059f6b74SHong Zhang 1321059f6b74SHong Zhang if (allrows[i]) { 1322059f6b74SHong Zhang row = row+rstart; 1323059f6b74SHong Zhang } else { 1324a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1325a42ab0d6SHong Zhang ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1326a42ab0d6SHong Zhang row--; 1327121971b2SHong Zhang 1328a42ab0d6SHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1329a42ab0d6SHong Zhang #else 1330a42ab0d6SHong Zhang row = rmap_i[row + rstart]; 1331a42ab0d6SHong Zhang #endif 1332059f6b74SHong Zhang } 1333a42ab0d6SHong Zhang mat_i = imat_i[row]; 1334a42ab0d6SHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1335a42ab0d6SHong Zhang mat_j = imat_j + mat_i; 133680d31651SHong Zhang ilen = imat_ilen[row]; 1337a42ab0d6SHong Zhang 1338a42ab0d6SHong Zhang /* load the column indices for this row into cols*/ 1339a42ab0d6SHong Zhang if (!allcolumns[i]) { 1340a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1341a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1342a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1343a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1344a42ab0d6SHong Zhang if (tcol) { 1345a42ab0d6SHong Zhang #else 1346a42ab0d6SHong Zhang if ((tcol = cmap_i[ctmp])) { 1347a42ab0d6SHong Zhang #endif 1348a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1349580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1350a42ab0d6SHong Zhang mat_a += bs2; 135180d31651SHong Zhang ilen++; 1352a42ab0d6SHong Zhang } 1353a42ab0d6SHong Zhang } else break; 1354a42ab0d6SHong Zhang } 1355a42ab0d6SHong Zhang imark = l; 1356a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1357a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1358a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1359a42ab0d6SHong Zhang if (tcol) { 1360a42ab0d6SHong Zhang #else 1361a42ab0d6SHong Zhang if ((tcol = cmap_i[cstart + cworkA[l]])) { 1362a42ab0d6SHong Zhang #endif 1363a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1364a42ab0d6SHong Zhang if (!ijonly) { 1365580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkA+l*bs2,bs2);CHKERRQ(ierr); 1366a42ab0d6SHong Zhang mat_a += bs2; 1367a42ab0d6SHong Zhang } 136880d31651SHong Zhang ilen++; 1369a42ab0d6SHong Zhang } 1370a42ab0d6SHong Zhang } 1371a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1372a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1373a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1374a42ab0d6SHong Zhang if (tcol) { 1375a42ab0d6SHong Zhang #else 1376a42ab0d6SHong Zhang if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1377a42ab0d6SHong Zhang #endif 1378a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1379a42ab0d6SHong Zhang if (!ijonly) { 1380580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1381a42ab0d6SHong Zhang mat_a += bs2; 1382a42ab0d6SHong Zhang } 138380d31651SHong Zhang ilen++; 1384a42ab0d6SHong Zhang } 1385a42ab0d6SHong Zhang } 1386a42ab0d6SHong Zhang } else { /* allcolumns */ 1387a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1388a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1389a42ab0d6SHong Zhang *mat_j++ = ctmp; 1390580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1391a42ab0d6SHong Zhang mat_a += bs2; 139280d31651SHong Zhang ilen++; 1393a42ab0d6SHong Zhang } else break; 1394a42ab0d6SHong Zhang } 1395a42ab0d6SHong Zhang imark = l; 1396a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1397a42ab0d6SHong Zhang *mat_j++ = cstart+cworkA[l]; 1398a42ab0d6SHong Zhang if (!ijonly) { 1399580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkA+l*bs2,bs2);CHKERRQ(ierr); 1400a42ab0d6SHong Zhang mat_a += bs2; 1401a42ab0d6SHong Zhang } 140280d31651SHong Zhang ilen++; 1403a42ab0d6SHong Zhang } 1404a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1405a42ab0d6SHong Zhang *mat_j++ = bmap[cworkB[l]]; 1406a42ab0d6SHong Zhang if (!ijonly) { 1407580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1408a42ab0d6SHong Zhang mat_a += bs2; 1409a42ab0d6SHong Zhang } 141080d31651SHong Zhang ilen++; 1411a42ab0d6SHong Zhang } 1412a42ab0d6SHong Zhang } 141380d31651SHong Zhang imat_ilen[row] = ilen; 141417df9f7cSHong Zhang } 141517df9f7cSHong Zhang } 141617df9f7cSHong Zhang } 141717df9f7cSHong Zhang 141817df9f7cSHong Zhang /* Now assemble the off proc rows */ 14195dd0c08cSHong Zhang if (!ijonly) { 142017df9f7cSHong Zhang ierr = MPI_Waitall(nrqs,r_waits4,r_status4);CHKERRQ(ierr); 14215dd0c08cSHong Zhang } 142217df9f7cSHong Zhang for (tmp2=0; tmp2<nrqs; tmp2++) { 142317df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 142417df9f7cSHong Zhang jmax = sbuf1_i[0]; 142517df9f7cSHong Zhang ct1 = 2*jmax + 1; 142617df9f7cSHong Zhang ct2 = 0; 142717df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 142817df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 14295dd0c08cSHong Zhang if (!ijonly) rbuf4_i = rbuf4[tmp2]; 143017df9f7cSHong Zhang for (j=1; j<=jmax; j++) { 143117df9f7cSHong Zhang is_no = sbuf1_i[2*j-1]; 143217df9f7cSHong Zhang rmap_i = rmap[is_no]; 143317df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 143417df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[is_no]->data; 143517df9f7cSHong Zhang imat_ilen = subc->ilen; 143617df9f7cSHong Zhang imat_j = subc->j; 143717df9f7cSHong Zhang imat_i = subc->i; 14385dd0c08cSHong Zhang if (!ijonly) imat_a = subc->a; 143917df9f7cSHong Zhang max1 = sbuf1_i[2*j]; 14405dd0c08cSHong Zhang for (k=0; k<max1; k++,ct1++) { /* for each recved block row */ 144117df9f7cSHong Zhang row = sbuf1_i[ct1]; 1442059f6b74SHong Zhang 1443059f6b74SHong Zhang if (allrows[is_no]) { 1444059f6b74SHong Zhang row = sbuf1_i[ct1]; 1445059f6b74SHong Zhang } else { 144617df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 144717df9f7cSHong Zhang ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 144817df9f7cSHong Zhang row--; 14495dd0c08cSHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 145017df9f7cSHong Zhang #else 145117df9f7cSHong Zhang row = rmap_i[row]; 145217df9f7cSHong Zhang #endif 1453059f6b74SHong Zhang } 145417df9f7cSHong Zhang ilen = imat_ilen[row]; 145517df9f7cSHong Zhang mat_i = imat_i[row]; 14565dd0c08cSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 145717df9f7cSHong Zhang mat_j = imat_j + mat_i; 145817df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 145917df9f7cSHong Zhang if (!allcolumns[is_no]) { 146017df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 146117df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 146217df9f7cSHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 146317df9f7cSHong Zhang #else 146417df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 146517df9f7cSHong Zhang #endif 146617df9f7cSHong Zhang if (tcol) { 146717df9f7cSHong Zhang *mat_j++ = tcol - 1; 14685dd0c08cSHong Zhang if (!ijonly) { 1469580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,rbuf4_i+ct2*bs2,bs2);CHKERRQ(ierr); 14705dd0c08cSHong Zhang mat_a += bs2; 14715dd0c08cSHong Zhang } 147217df9f7cSHong Zhang ilen++; 147317df9f7cSHong Zhang } 147417df9f7cSHong Zhang } 147517df9f7cSHong Zhang } else { /* allcolumns */ 147617df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 147717df9f7cSHong Zhang *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 14785dd0c08cSHong Zhang if (!ijonly) { 1479580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,rbuf4_i+ct2*bs2,bs2);CHKERRQ(ierr); 14805dd0c08cSHong Zhang mat_a += bs2; 14815dd0c08cSHong Zhang } 148217df9f7cSHong Zhang ilen++; 148317df9f7cSHong Zhang } 148417df9f7cSHong Zhang } 148517df9f7cSHong Zhang imat_ilen[row] = ilen; 148617df9f7cSHong Zhang } 148717df9f7cSHong Zhang } 148817df9f7cSHong Zhang } 148917df9f7cSHong Zhang 149080d31651SHong Zhang if (!iscsorted) { /* sort column indices of the rows */ 149180d31651SHong Zhang MatScalar *work; 149280d31651SHong Zhang 149380d31651SHong Zhang ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 149417df9f7cSHong Zhang for (i=0; i<ismax; i++) { 149517df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 149680d31651SHong Zhang imat_ilen = subc->ilen; 149717df9f7cSHong Zhang imat_j = subc->j; 149817df9f7cSHong Zhang imat_i = subc->i; 14994b6ceb0dSHong Zhang if (!ijonly) imat_a = subc->a; 150017df9f7cSHong Zhang if (allcolumns[i]) continue; 15014b6ceb0dSHong Zhang 150217df9f7cSHong Zhang jmax = nrow[i]; 150317df9f7cSHong Zhang for (j=0; j<jmax; j++) { 150480d31651SHong Zhang mat_i = imat_i[j]; 150580d31651SHong Zhang mat_j = imat_j + mat_i; 15064b6ceb0dSHong Zhang ilen = imat_ilen[j]; 150780d31651SHong Zhang if (ijonly) { 150880d31651SHong Zhang ierr = PetscSortInt(ilen,mat_j);CHKERRQ(ierr); 150980d31651SHong Zhang } else { 15104b6ceb0dSHong Zhang mat_a = imat_a + mat_i*bs2; 151180d31651SHong Zhang ierr = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr); 151217df9f7cSHong Zhang } 151317df9f7cSHong Zhang } 151417df9f7cSHong Zhang } 151580d31651SHong Zhang ierr = PetscFree(work);CHKERRQ(ierr); 151680d31651SHong Zhang } 151717df9f7cSHong Zhang 1518bbc89d27SHong Zhang if (!ijonly) { 151917df9f7cSHong Zhang ierr = PetscFree(r_status4);CHKERRQ(ierr); 152017df9f7cSHong Zhang ierr = PetscFree(r_waits4);CHKERRQ(ierr); 152117df9f7cSHong Zhang if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 152217df9f7cSHong Zhang ierr = PetscFree(s_waits4);CHKERRQ(ierr); 152317df9f7cSHong Zhang ierr = PetscFree(s_status4);CHKERRQ(ierr); 1524bbc89d27SHong Zhang } 152517df9f7cSHong Zhang 152617df9f7cSHong Zhang /* Restore the indices */ 152717df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1528ec3c739cSHong Zhang if (!allrows[i]) { 152917df9f7cSHong Zhang ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1530ec3c739cSHong Zhang } 153117df9f7cSHong Zhang if (!allcolumns[i]) { 153217df9f7cSHong Zhang ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 153317df9f7cSHong Zhang } 153417df9f7cSHong Zhang } 153517df9f7cSHong Zhang 153617df9f7cSHong Zhang for (i=0; i<ismax; i++) { 153717df9f7cSHong Zhang ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153817df9f7cSHong Zhang ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153917df9f7cSHong Zhang } 154017df9f7cSHong Zhang 15412cd485c2SBarry Smith ierr = PetscFree5(*(PetscInt***)&irow,*(PetscInt***)&icol,nrow,ncol,issorted);CHKERRQ(ierr); 15424698041cSHong Zhang ierr = PetscFree5(row2proc,cmap,rmap,allcolumns,allrows);CHKERRQ(ierr); 154317df9f7cSHong Zhang 1544bbc89d27SHong Zhang if (!ijonly) { 154517df9f7cSHong Zhang ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 154617df9f7cSHong Zhang ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 154717df9f7cSHong Zhang 154817df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 154917df9f7cSHong Zhang ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 155017df9f7cSHong Zhang } 155117df9f7cSHong Zhang ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1552bbc89d27SHong Zhang } 15537a868f3eSHong Zhang c->ijonly = PETSC_FALSE; /* set back to the default */ 15543a40ed3dSBarry Smith PetscFunctionReturn(0); 1555d5b485abSSatish Balay } 1556