1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 30500aeceSBarry Smith /* 40500aeceSBarry Smith Support for the parallel dense matrix vector multiply 50500aeceSBarry Smith */ 670f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 70500aeceSBarry Smith 84a2ae208SSatish Balay #undef __FUNCT__ 94a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIDense" 10dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIDense(Mat mat) 110500aeceSBarry Smith { 120500aeceSBarry Smith Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 13dfbe8321SBarry Smith PetscErrorCode ierr; 14b9b97703SBarry Smith IS from,to; 150500aeceSBarry Smith Vec gvec; 160500aeceSBarry Smith 173a40ed3dSBarry Smith PetscFunctionBegin; 180500aeceSBarry Smith /* Create local vector that is used to scatter into */ 19d0f46423SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,mat->cmap->N,&mdn->lvec);CHKERRQ(ierr); 200500aeceSBarry Smith 210500aeceSBarry Smith /* Create temporary index set for building scatter gather */ 22d0f46423SBarry Smith ierr = ISCreateStride(((PetscObject)mat)->comm,mat->cmap->N,0,1,&from);CHKERRQ(ierr); 23d0f46423SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,mat->cmap->N,0,1,&to);CHKERRQ(ierr); 240500aeceSBarry Smith 250500aeceSBarry Smith /* Create temporary global vector to generate scatter context */ 266f0b7f37SSatish Balay /* n = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank]; */ 27426b7fd2SBarry Smith 28d0f46423SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,mdn->nvec,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr); 29e19c7942SLois Curfman McInnes 300500aeceSBarry Smith /* Generate the scatter context */ 31b9b97703SBarry Smith ierr = VecScatterCreate(gvec,from,mdn->lvec,to,&mdn->Mvctx);CHKERRQ(ierr); 3252e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,mdn->Mvctx);CHKERRQ(ierr); 3352e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,mdn->lvec);CHKERRQ(ierr); 3452e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 3552e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 3652e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,gvec);CHKERRQ(ierr); 37e19c7942SLois Curfman McInnes 38b9b97703SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 39b9b97703SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 40d9f96c7cSLois Curfman McInnes ierr = VecDestroy(gvec);CHKERRQ(ierr); 413a40ed3dSBarry Smith PetscFunctionReturn(0); 420500aeceSBarry Smith } 430500aeceSBarry Smith 4413f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_MPIDense_Local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*); 454a2ae208SSatish Balay #undef __FUNCT__ 464a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIDense" 4713f74950SBarry Smith PetscErrorCode MatGetSubMatrices_MPIDense(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 48e26ad7d8SSatish Balay { 496849ba73SBarry Smith PetscErrorCode ierr; 5013f74950SBarry Smith PetscInt nmax,nstages_local,nstages,i,pos,max_no; 51e26ad7d8SSatish Balay 52e26ad7d8SSatish Balay PetscFunctionBegin; 53e26ad7d8SSatish Balay /* Allocate memory to hold all the submatrices */ 54e26ad7d8SSatish Balay if (scall != MAT_REUSE_MATRIX) { 55b0a32e0cSBarry Smith ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr); 56e26ad7d8SSatish Balay } 57e26ad7d8SSatish Balay /* Determine the number of stages through which submatrices are done */ 58d0f46423SBarry Smith nmax = 20*1000000 / (C->cmap->N * sizeof(PetscInt)); 59329f5518SBarry Smith if (!nmax) nmax = 1; 60e26ad7d8SSatish Balay nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 61e26ad7d8SSatish Balay 62e26ad7d8SSatish Balay /* Make sure every processor loops through the nstages */ 637adad957SLisandro Dalcin ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);CHKERRQ(ierr); 640500aeceSBarry Smith 650500aeceSBarry Smith 66e26ad7d8SSatish Balay for (i=0,pos=0; i<nstages; i++) { 67e26ad7d8SSatish Balay if (pos+nmax <= ismax) max_no = nmax; 68e26ad7d8SSatish Balay else if (pos == ismax) max_no = 0; 69e26ad7d8SSatish Balay else max_no = ismax-pos; 70e26ad7d8SSatish Balay ierr = MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);CHKERRQ(ierr); 71e26ad7d8SSatish Balay pos += max_no; 72e26ad7d8SSatish Balay } 73e26ad7d8SSatish Balay PetscFunctionReturn(0); 74e26ad7d8SSatish Balay } 75e26ad7d8SSatish Balay /* -------------------------------------------------------------------------*/ 764a2ae208SSatish Balay #undef __FUNCT__ 774a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIDense_Local" 7813f74950SBarry Smith PetscErrorCode MatGetSubMatrices_MPIDense_Local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats) 79e26ad7d8SSatish Balay { 80e26ad7d8SSatish Balay Mat_MPIDense *c = (Mat_MPIDense*)C->data; 81e26ad7d8SSatish Balay Mat A = c->A; 82e26ad7d8SSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data,*mat; 836849ba73SBarry Smith PetscErrorCode ierr; 8413f74950SBarry Smith PetscMPIInt rank,size,tag0,tag1,idex,end,i; 85d0f46423SBarry Smith PetscInt N = C->cmap->N,rstart = C->rmap->rstart,count; 86*5d0c19d7SBarry Smith const PetscInt **irow,**icol,*irow_i; 87*5d0c19d7SBarry Smith PetscInt *nrow,*ncol,*w1,*w3,*w4,*rtable,start; 8813f74950SBarry Smith PetscInt **sbuf1,m,j,k,l,ct1,**rbuf1,row,proc; 8913f74950SBarry Smith PetscInt nrqs,msz,**ptr,*ctr,*pa,*tmp,bsz,nrqr; 90*5d0c19d7SBarry Smith PetscInt is_no,jmax,**rmap,*rmap_i; 9113f74950SBarry Smith PetscInt len,ctr_j,*sbuf1_j,*rbuf1_i; 924b3ec233SSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 934b3ec233SSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status2; 94e26ad7d8SSatish Balay MPI_Comm comm; 9587828ca2SBarry Smith PetscScalar **rbuf2,**sbuf2; 965c5985e7SKris Buschelman PetscTruth sorted; 97e26ad7d8SSatish Balay 98e26ad7d8SSatish Balay PetscFunctionBegin; 997adad957SLisandro Dalcin comm = ((PetscObject)C)->comm; 1007adad957SLisandro Dalcin tag0 = ((PetscObject)C)->tag; 101e26ad7d8SSatish Balay size = c->size; 102e26ad7d8SSatish Balay rank = c->rank; 103d0f46423SBarry Smith m = C->rmap->N; 104e26ad7d8SSatish Balay 105e26ad7d8SSatish Balay /* Get some new tags to keep the communication clean */ 106e26ad7d8SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 107e26ad7d8SSatish Balay 108e26ad7d8SSatish Balay /* Check if the col indices are sorted */ 109e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 1105c5985e7SKris Buschelman ierr = ISSorted(isrow[i],&sorted);CHKERRQ(ierr); 1115c5985e7SKris Buschelman if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted"); 1125c5985e7SKris Buschelman ierr = ISSorted(iscol[i],&sorted);CHKERRQ(ierr); 1135c5985e7SKris Buschelman if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted"); 114e26ad7d8SSatish Balay } 115e26ad7d8SSatish Balay 11613f74950SBarry Smith len = 2*ismax*(sizeof(PetscInt*)+sizeof(PetscInt)) + (m+1)*sizeof(PetscInt); 117b0a32e0cSBarry Smith ierr = PetscMalloc(len,&irow);CHKERRQ(ierr); 118e26ad7d8SSatish Balay icol = irow + ismax; 11913f74950SBarry Smith nrow = (PetscInt*)(icol + ismax); 120e26ad7d8SSatish Balay ncol = nrow + ismax; 121e26ad7d8SSatish Balay rtable = ncol + ismax; 122e26ad7d8SSatish Balay 123e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 124e26ad7d8SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 125e26ad7d8SSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 126b9b97703SBarry Smith ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 127b9b97703SBarry Smith ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 128e26ad7d8SSatish Balay } 129e26ad7d8SSatish Balay 130e26ad7d8SSatish Balay /* Create hash table for the mapping :row -> proc*/ 131e26ad7d8SSatish Balay for (i=0,j=0; i<size; i++) { 132d0f46423SBarry Smith jmax = C->rmap->range[i+1]; 133e26ad7d8SSatish Balay for (; j<jmax; j++) { 134e26ad7d8SSatish Balay rtable[j] = i; 135e26ad7d8SSatish Balay } 136e26ad7d8SSatish Balay } 137e26ad7d8SSatish Balay 138e26ad7d8SSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 139e26ad7d8SSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 14013f74950SBarry Smith ierr = PetscMalloc(size*4*sizeof(PetscInt),&w1);CHKERRQ(ierr); /* mesg size */ 141c1dc657dSBarry Smith w3 = w1 + 2*size; /* no of IS that needs to be sent to proc i */ 142c1dc657dSBarry Smith w4 = w3 + size; /* temp work space used in determining w1, w3 */ 14313f74950SBarry Smith ierr = PetscMemzero(w1,size*3*sizeof(PetscInt));CHKERRQ(ierr); /* initialize work vector*/ 144e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 14513f74950SBarry Smith ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialize work vector*/ 146e26ad7d8SSatish Balay jmax = nrow[i]; 147e26ad7d8SSatish Balay irow_i = irow[i]; 148e26ad7d8SSatish Balay for (j=0; j<jmax; j++) { 149e26ad7d8SSatish Balay row = irow_i[j]; 150e26ad7d8SSatish Balay proc = rtable[row]; 151e26ad7d8SSatish Balay w4[proc]++; 152e26ad7d8SSatish Balay } 153e26ad7d8SSatish Balay for (j=0; j<size; j++) { 154c1dc657dSBarry Smith if (w4[j]) { w1[2*j] += w4[j]; w3[j]++;} 155e26ad7d8SSatish Balay } 156e26ad7d8SSatish Balay } 157e26ad7d8SSatish Balay 158e26ad7d8SSatish Balay nrqs = 0; /* no of outgoing messages */ 159e26ad7d8SSatish Balay msz = 0; /* total mesg length (for all procs) */ 160c1dc657dSBarry Smith w1[2*rank] = 0; /* no mesg sent to self */ 161e26ad7d8SSatish Balay w3[rank] = 0; 162e26ad7d8SSatish Balay for (i=0; i<size; i++) { 163c1dc657dSBarry Smith if (w1[2*i]) { w1[2*i+1] = 1; nrqs++;} /* there exists a message to proc i */ 164e26ad7d8SSatish Balay } 16513f74950SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/ 166e26ad7d8SSatish Balay for (i=0,j=0; i<size; i++) { 167c1dc657dSBarry Smith if (w1[2*i]) { pa[j] = i; j++; } 168e26ad7d8SSatish Balay } 169e26ad7d8SSatish Balay 170e26ad7d8SSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 171e26ad7d8SSatish Balay for (i=0; i<nrqs; i++) { 172e26ad7d8SSatish Balay j = pa[i]; 173c1dc657dSBarry Smith w1[2*j] += w1[2*j+1] + 2* w3[j]; 174c1dc657dSBarry Smith msz += w1[2*j]; 175e26ad7d8SSatish Balay } 176e26ad7d8SSatish Balay /* Do a global reduction to determine how many messages to expect*/ 177c1dc657dSBarry Smith ierr = PetscMaxSum(comm,w1,&bsz,&nrqr);CHKERRQ(ierr); 178e26ad7d8SSatish Balay 179e26ad7d8SSatish Balay /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 18013f74950SBarry Smith len = (nrqr+1)*sizeof(PetscInt*) + nrqr*bsz*sizeof(PetscInt); 181b0a32e0cSBarry Smith ierr = PetscMalloc(len,&rbuf1);CHKERRQ(ierr); 18213f74950SBarry Smith rbuf1[0] = (PetscInt*)(rbuf1 + nrqr); 183e26ad7d8SSatish Balay for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz; 184e26ad7d8SSatish Balay 185e26ad7d8SSatish Balay /* Post the receives */ 186b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&r_waits1);CHKERRQ(ierr); 187e26ad7d8SSatish Balay for (i=0; i<nrqr; ++i) { 18813f74950SBarry Smith ierr = MPI_Irecv(rbuf1[i],bsz,MPIU_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr); 189e26ad7d8SSatish Balay } 190e26ad7d8SSatish Balay 191e26ad7d8SSatish Balay /* Allocate Memory for outgoing messages */ 19213f74950SBarry Smith len = 2*size*sizeof(PetscInt*) + 2*msz*sizeof(PetscInt)+ size*sizeof(PetscInt); 193b0a32e0cSBarry Smith ierr = PetscMalloc(len,&sbuf1);CHKERRQ(ierr); 194e26ad7d8SSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 19513f74950SBarry Smith ierr = PetscMemzero(sbuf1,2*size*sizeof(PetscInt*));CHKERRQ(ierr); 196e26ad7d8SSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 19713f74950SBarry Smith tmp = (PetscInt*)(ptr + size); 198e26ad7d8SSatish Balay ctr = tmp + 2*msz; 199e26ad7d8SSatish Balay 200e26ad7d8SSatish Balay { 20113f74950SBarry Smith PetscInt *iptr = tmp,ict = 0; 202e26ad7d8SSatish Balay for (i=0; i<nrqs; i++) { 203e26ad7d8SSatish Balay j = pa[i]; 204e26ad7d8SSatish Balay iptr += ict; 205e26ad7d8SSatish Balay sbuf1[j] = iptr; 206c1dc657dSBarry Smith ict = w1[2*j]; 207e26ad7d8SSatish Balay } 208e26ad7d8SSatish Balay } 209e26ad7d8SSatish Balay 210e26ad7d8SSatish Balay /* Form the outgoing messages */ 211e26ad7d8SSatish Balay /* Initialize the header space */ 212e26ad7d8SSatish Balay for (i=0; i<nrqs; i++) { 213e26ad7d8SSatish Balay j = pa[i]; 214e26ad7d8SSatish Balay sbuf1[j][0] = 0; 21513f74950SBarry Smith ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 216e26ad7d8SSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 217e26ad7d8SSatish Balay } 218e26ad7d8SSatish Balay 219e26ad7d8SSatish Balay /* Parse the isrow and copy data into outbuf */ 220e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 22113f74950SBarry Smith ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 222e26ad7d8SSatish Balay irow_i = irow[i]; 223e26ad7d8SSatish Balay jmax = nrow[i]; 224e26ad7d8SSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 225e26ad7d8SSatish Balay row = irow_i[j]; 226e26ad7d8SSatish Balay proc = rtable[row]; 227e26ad7d8SSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 228e26ad7d8SSatish Balay ctr[proc]++; 229e26ad7d8SSatish Balay *ptr[proc] = row; 230e26ad7d8SSatish Balay ptr[proc]++; 231e26ad7d8SSatish Balay } 232e26ad7d8SSatish Balay } 233e26ad7d8SSatish Balay /* Update the headers for the current IS */ 234e26ad7d8SSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 23506ef35abSSatish Balay if ((ctr_j = ctr[j])) { 236e26ad7d8SSatish Balay sbuf1_j = sbuf1[j]; 237e26ad7d8SSatish Balay k = ++sbuf1_j[0]; 238e26ad7d8SSatish Balay sbuf1_j[2*k] = ctr_j; 239e26ad7d8SSatish Balay sbuf1_j[2*k-1] = i; 240e26ad7d8SSatish Balay } 241e26ad7d8SSatish Balay } 242e26ad7d8SSatish Balay } 243e26ad7d8SSatish Balay 244e26ad7d8SSatish Balay /* Now post the sends */ 245b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 246e26ad7d8SSatish Balay for (i=0; i<nrqs; ++i) { 247e26ad7d8SSatish Balay j = pa[i]; 24813f74950SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[2*j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 249e26ad7d8SSatish Balay } 250e26ad7d8SSatish Balay 2514b3ec233SSatish Balay /* Post recieves to capture the row_data from other procs */ 252b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr); 25387828ca2SBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(PetscScalar*),&rbuf2);CHKERRQ(ierr); 2544b3ec233SSatish Balay for (i=0; i<nrqs; i++) { 255e26ad7d8SSatish Balay j = pa[i]; 256c1dc657dSBarry Smith count = (w1[2*j] - (2*sbuf1[j][0] + 1))*N; 25787828ca2SBarry Smith ierr = PetscMalloc((count+1)*sizeof(PetscScalar),&rbuf2[i]);CHKERRQ(ierr); 2584b3ec233SSatish Balay ierr = MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 259e26ad7d8SSatish Balay } 260e26ad7d8SSatish Balay 2614b3ec233SSatish Balay /* Receive messages(row_nos) and then, pack and send off the rowvalues 2624b3ec233SSatish Balay to the correct processors */ 263e26ad7d8SSatish Balay 264b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 265b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr); 26687828ca2SBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(PetscScalar*),&sbuf2);CHKERRQ(ierr); 267e26ad7d8SSatish Balay 268e26ad7d8SSatish Balay { 26987828ca2SBarry Smith PetscScalar *sbuf2_i,*v_start; 27013f74950SBarry Smith PetscInt s_proc; 271e26ad7d8SSatish Balay for (i=0; i<nrqr; ++i) { 272999d9058SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 2734b3ec233SSatish Balay s_proc = r_status1[i].MPI_SOURCE; /* send processor */ 274999d9058SBarry Smith rbuf1_i = rbuf1[idex]; /* Actual message from s_proc */ 275999d9058SBarry Smith /* no of rows = end - start; since start is array idex[], 0idex, whel end 276999d9058SBarry Smith is length of the buffer - which is 1idex */ 277e26ad7d8SSatish Balay start = 2*rbuf1_i[0] + 1; 27813f74950SBarry Smith ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 2794b3ec233SSatish Balay /* allocate memory sufficinet to hold all the row values */ 280999d9058SBarry Smith ierr = PetscMalloc((end-start)*N*sizeof(PetscScalar),&sbuf2[idex]);CHKERRQ(ierr); 281999d9058SBarry Smith sbuf2_i = sbuf2[idex]; 2824b3ec233SSatish Balay /* Now pack the data */ 283e26ad7d8SSatish Balay for (j=start; j<end; j++) { 2844b3ec233SSatish Balay row = rbuf1_i[j] - rstart; 2854b3ec233SSatish Balay v_start = a->v + row; 2864b3ec233SSatish Balay for (k=0; k<N; k++) { 2874b3ec233SSatish Balay sbuf2_i[0] = v_start[0]; 288d0f46423SBarry Smith sbuf2_i++; v_start += C->rmap->n; 289e26ad7d8SSatish Balay } 290e26ad7d8SSatish Balay } 2914b3ec233SSatish Balay /* Now send off the data */ 292999d9058SBarry Smith ierr = MPI_Isend(sbuf2[idex],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);CHKERRQ(ierr); 2934b3ec233SSatish Balay } 2944b3ec233SSatish Balay } 2954b3ec233SSatish Balay /* End Send-Recv of IS + row_numbers */ 296606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 297606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 298b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr); 2990c468ba9SBarry Smith if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 300606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 301606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 302e26ad7d8SSatish Balay 3034b3ec233SSatish Balay /* Create the submatrices */ 3044b3ec233SSatish Balay if (scall == MAT_REUSE_MATRIX) { 305e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 3064b3ec233SSatish Balay mat = (Mat_SeqDense *)(submats[i]->data); 307d0f46423SBarry Smith if ((submats[i]->rmap->n != nrow[i]) || (submats[i]->cmap->n != ncol[i])) { 30829bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 309e26ad7d8SSatish Balay } 310d0f46423SBarry Smith ierr = PetscMemzero(mat->v,submats[i]->rmap->n*submats[i]->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 3114b3ec233SSatish Balay submats[i]->factor = C->factor; 312e26ad7d8SSatish Balay } 3134b3ec233SSatish Balay } else { 314e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 315f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 316f69a0ea3SMatthew Knepley ierr = MatSetSizes(submats[i],nrow[i],ncol[i],nrow[i],ncol[i]);CHKERRQ(ierr); 3177adad957SLisandro Dalcin ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 3185c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(submats[i],PETSC_NULL);CHKERRQ(ierr); 3194b3ec233SSatish Balay } 3204b3ec233SSatish Balay } 3214b3ec233SSatish Balay 3224b3ec233SSatish Balay /* Assemble the matrices */ 3234b3ec233SSatish Balay { 32413f74950SBarry Smith PetscInt col; 32587828ca2SBarry Smith PetscScalar *imat_v,*mat_v,*imat_vi,*mat_vi; 3264b3ec233SSatish Balay 3274b3ec233SSatish Balay for (i=0; i<ismax; i++) { 3284b3ec233SSatish Balay mat = (Mat_SeqDense*)submats[i]->data; 3294b3ec233SSatish Balay mat_v = a->v; 3304b3ec233SSatish Balay imat_v = mat->v; 331e26ad7d8SSatish Balay irow_i = irow[i]; 3324b3ec233SSatish Balay m = nrow[i]; 3334b3ec233SSatish Balay for (j=0; j<m; j++) { 334e26ad7d8SSatish Balay row = irow_i[j] ; 335e26ad7d8SSatish Balay proc = rtable[row]; 336e26ad7d8SSatish Balay if (proc == rank) { 3374b3ec233SSatish Balay row = row - rstart; 3384b3ec233SSatish Balay mat_vi = mat_v + row; 3394b3ec233SSatish Balay imat_vi = imat_v + j; 3404b3ec233SSatish Balay for (k=0; k<ncol[i]; k++) { 3414b3ec233SSatish Balay col = icol[i][k]; 342d0f46423SBarry Smith imat_vi[k*m] = mat_vi[col*C->rmap->n]; 343e26ad7d8SSatish Balay } 3444b3ec233SSatish Balay } 345e26ad7d8SSatish Balay } 346e26ad7d8SSatish Balay } 347e26ad7d8SSatish Balay } 348e26ad7d8SSatish Balay 349d0f46423SBarry Smith /* Create row map-> This maps c->row to submat->row for each submat*/ 3504b3ec233SSatish Balay /* this is a very expensive operation wrt memory usage */ 351d0f46423SBarry Smith len = (1+ismax)*sizeof(PetscInt*)+ ismax*C->rmap->N*sizeof(PetscInt); 352b0a32e0cSBarry Smith ierr = PetscMalloc(len,&rmap);CHKERRQ(ierr); 35313f74950SBarry Smith rmap[0] = (PetscInt*)(rmap + ismax); 354d0f46423SBarry Smith ierr = PetscMemzero(rmap[0],ismax*C->rmap->N*sizeof(PetscInt));CHKERRQ(ierr); 355d0f46423SBarry Smith for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->rmap->N;} 356e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 357e26ad7d8SSatish Balay rmap_i = rmap[i]; 358e26ad7d8SSatish Balay irow_i = irow[i]; 359e26ad7d8SSatish Balay jmax = nrow[i]; 360e26ad7d8SSatish Balay for (j=0; j<jmax; j++) { 361e26ad7d8SSatish Balay rmap_i[irow_i[j]] = j; 362e26ad7d8SSatish Balay } 363e26ad7d8SSatish Balay } 364e26ad7d8SSatish Balay 3654b3ec233SSatish Balay /* Now Receive the row_values and assemble the rest of the matrix */ 366b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr); 3674b3ec233SSatish Balay 3684b3ec233SSatish Balay { 36913f74950SBarry Smith PetscInt is_max,tmp1,col,*sbuf1_i,is_sz; 37087828ca2SBarry Smith PetscScalar *rbuf2_i,*imat_v,*imat_vi; 3714b3ec233SSatish Balay 3724b3ec233SSatish Balay for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */ 3734b3ec233SSatish Balay ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr); 3744b3ec233SSatish Balay /* Now dig out the corresponding sbuf1, which contains the IS data_structure */ 3754b3ec233SSatish Balay sbuf1_i = sbuf1[pa[i]]; 3764b3ec233SSatish Balay is_max = sbuf1_i[0]; 3774b3ec233SSatish Balay ct1 = 2*is_max+1; 378e26ad7d8SSatish Balay rbuf2_i = rbuf2[i]; 3794b3ec233SSatish Balay for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */ 380e26ad7d8SSatish Balay is_no = sbuf1_i[2*j-1]; 3814b3ec233SSatish Balay is_sz = sbuf1_i[2*j]; 382e26ad7d8SSatish Balay mat = (Mat_SeqDense*)submats[is_no]->data; 3834b3ec233SSatish Balay imat_v = mat->v; 3844b3ec233SSatish Balay rmap_i = rmap[is_no]; 3854b3ec233SSatish Balay m = nrow[is_no]; 3864b3ec233SSatish Balay for (k=0; k<is_sz; k++,rbuf2_i+=N) { /* For each row */ 3874b3ec233SSatish Balay row = sbuf1_i[ct1]; ct1++; 388e26ad7d8SSatish Balay row = rmap_i[row]; 3894b3ec233SSatish Balay imat_vi = imat_v + row; 3904b3ec233SSatish Balay for (l=0; l<ncol[is_no]; l++) { /* For each col */ 3914b3ec233SSatish Balay col = icol[is_no][l]; 3924b3ec233SSatish Balay imat_vi[l*m] = rbuf2_i[col]; 393e26ad7d8SSatish Balay } 394e26ad7d8SSatish Balay } 395e26ad7d8SSatish Balay } 396e26ad7d8SSatish Balay } 3974b3ec233SSatish Balay } 3984b3ec233SSatish Balay /* End Send-Recv of row_values */ 399606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 400606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 401b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr); 4020c468ba9SBarry Smith if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 403606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 404606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 405e26ad7d8SSatish Balay 406e26ad7d8SSatish Balay /* Restore the indices */ 407e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 408e26ad7d8SSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 409e26ad7d8SSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 410e26ad7d8SSatish Balay } 411e26ad7d8SSatish Balay 412e26ad7d8SSatish Balay /* Destroy allocated memory */ 413606d414cSSatish Balay ierr = PetscFree(irow);CHKERRQ(ierr); 414606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 415606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 416e26ad7d8SSatish Balay 4174b3ec233SSatish Balay 4184b3ec233SSatish Balay for (i=0; i<nrqs; ++i) { 419606d414cSSatish Balay ierr = PetscFree(rbuf2[i]);CHKERRQ(ierr); 4204b3ec233SSatish Balay } 421606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 422606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 423606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 4244b3ec233SSatish Balay 425e26ad7d8SSatish Balay for (i=0; i<nrqr; ++i) { 426606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 427e26ad7d8SSatish Balay } 428e26ad7d8SSatish Balay 429606d414cSSatish Balay ierr = PetscFree(sbuf2);CHKERRQ(ierr); 430606d414cSSatish Balay ierr = PetscFree(rmap);CHKERRQ(ierr); 431e26ad7d8SSatish Balay 432e26ad7d8SSatish Balay for (i=0; i<ismax; i++) { 433e26ad7d8SSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 434e26ad7d8SSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 435e26ad7d8SSatish Balay } 436e26ad7d8SSatish Balay 437e26ad7d8SSatish Balay PetscFunctionReturn(0); 438e26ad7d8SSatish Balay } 439e26ad7d8SSatish Balay 440