1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*be0abb6dSBarry Smith static char vcid[] = "$Id: mmdense.c,v 1.23 1999/08/17 19:00:49 balay Exp bsmith $"; 30500aeceSBarry Smith #endif 40500aeceSBarry Smith 50500aeceSBarry Smith /* 60500aeceSBarry Smith Support for the parallel dense matrix vector multiply 70500aeceSBarry Smith */ 870f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 9f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 100500aeceSBarry Smith 115615d1e5SSatish Balay #undef __FUNC__ 125615d1e5SSatish Balay #define __FUNC__ "MatSetUpMultiply_MPIDense" 130500aeceSBarry Smith int MatSetUpMultiply_MPIDense(Mat mat) 140500aeceSBarry Smith { 150500aeceSBarry Smith Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 166f0b7f37SSatish Balay int ierr; 17d9f96c7cSLois Curfman McInnes IS tofrom; 180500aeceSBarry Smith Vec gvec; 190500aeceSBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 210500aeceSBarry Smith /* Create local vector that is used to scatter into */ 22029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,mdn->N,&mdn->lvec);CHKERRQ(ierr); 230500aeceSBarry Smith 240500aeceSBarry Smith /* Create temporary index set for building scatter gather */ 25029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,mdn->N,0,1,&tofrom);CHKERRQ(ierr); 260500aeceSBarry Smith 270500aeceSBarry Smith /* Create temporary global vector to generate scatter context */ 286f0b7f37SSatish Balay /* n = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank]; */ 29426b7fd2SBarry Smith 30*be0abb6dSBarry Smith ierr = VecCreateMPI(mat->comm,mdn->nvec,mdn->N,&gvec);CHKERRQ(ierr); 31e19c7942SLois Curfman McInnes 320500aeceSBarry Smith /* Generate the scatter context */ 33d9f96c7cSLois Curfman McInnes ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx);CHKERRQ(ierr); 340500aeceSBarry Smith PLogObjectParent(mat,mdn->Mvctx); 350500aeceSBarry Smith PLogObjectParent(mat,mdn->lvec); 36d9f96c7cSLois Curfman McInnes PLogObjectParent(mat,tofrom); 37d9f96c7cSLois Curfman McInnes PLogObjectParent(mat,gvec); 38e19c7942SLois Curfman McInnes 39d9f96c7cSLois Curfman McInnes ierr = ISDestroy(tofrom);CHKERRQ(ierr); 40d9f96c7cSLois Curfman McInnes ierr = VecDestroy(gvec);CHKERRQ(ierr); 413a40ed3dSBarry Smith PetscFunctionReturn(0); 420500aeceSBarry Smith } 430500aeceSBarry Smith 447b2a1423SBarry Smith extern int MatGetSubMatrices_MPIDense_Local(Mat,int,IS*,IS*,MatReuse,Mat*); 45e26ad7d8SSatish Balay #undef __FUNC__ 46e26ad7d8SSatish Balay #define __FUNC__ "MatGetSubMatrices_MPIDense" 47e26ad7d8SSatish Balay int MatGetSubMatrices_MPIDense(Mat C,int ismax,IS *isrow,IS *iscol, 487b2a1423SBarry Smith MatReuse scall,Mat **submat) 49e26ad7d8SSatish Balay { 50e26ad7d8SSatish Balay Mat_MPIDense *c = (Mat_MPIDense *) C->data; 51e26ad7d8SSatish Balay int nmax,nstages_local,nstages,i,pos,max_no,ierr; 52e26ad7d8SSatish Balay 53e26ad7d8SSatish Balay PetscFunctionBegin; 54e26ad7d8SSatish Balay /* Allocate memory to hold all the submatrices */ 55e26ad7d8SSatish Balay if (scall != MAT_REUSE_MATRIX) { 56e26ad7d8SSatish Balay *submat = (Mat *)PetscMalloc((ismax+1)*sizeof(Mat));CHKPTRQ(*submat); 57e26ad7d8SSatish Balay } 58e26ad7d8SSatish Balay /* Determine the number of stages through which submatrices are done */ 59e26ad7d8SSatish Balay nmax = 20*1000000 / (c->N * sizeof(int)); 60e26ad7d8SSatish Balay if (nmax == 0) nmax = 1; 61e26ad7d8SSatish Balay nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 62e26ad7d8SSatish Balay 63e26ad7d8SSatish Balay /* Make sure every processor loops through the nstages */ 64e26ad7d8SSatish Balay ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr); 650500aeceSBarry Smith 660500aeceSBarry Smith 67e26ad7d8SSatish Balay for ( i=0,pos=0; i<nstages; i++ ) { 68e26ad7d8SSatish Balay if (pos+nmax <= ismax) max_no = nmax; 69e26ad7d8SSatish Balay else if (pos == ismax) max_no = 0; 70e26ad7d8SSatish Balay else max_no = ismax-pos; 71e26ad7d8SSatish Balay ierr = MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);CHKERRQ(ierr); 72e26ad7d8SSatish Balay pos += max_no; 73e26ad7d8SSatish Balay } 74e26ad7d8SSatish Balay PetscFunctionReturn(0); 75e26ad7d8SSatish Balay } 76e26ad7d8SSatish Balay /* -------------------------------------------------------------------------*/ 77e26ad7d8SSatish Balay #undef __FUNC__ 78e26ad7d8SSatish Balay #define __FUNC__ "MatGetSubMatrices_MPIDense_Local" 797b2a1423SBarry Smith int MatGetSubMatrices_MPIDense_Local(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat *submats) 80e26ad7d8SSatish Balay { 81e26ad7d8SSatish Balay Mat_MPIDense *c = (Mat_MPIDense *) C->data; 82e26ad7d8SSatish Balay Mat A = c->A; 83e26ad7d8SSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data, *mat; 844b3ec233SSatish Balay int N = c->N, rstart = c->rstart,count; 85e26ad7d8SSatish Balay int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size; 864b3ec233SSatish Balay int **sbuf1, rank, m,i,j,k,l,ct1,ierr, **rbuf1,row,proc; 874b3ec233SSatish Balay int nrqs, msz, **ptr,index,*ctr,*pa,*tmp,bsz,nrqr; 884b3ec233SSatish Balay int is_no,jmax,*irow_i,**rmap,*rmap_i; 894b3ec233SSatish Balay int len,ctr_j,*sbuf1_j,*rbuf1_i; 904b3ec233SSatish Balay int tag0,tag1; 914b3ec233SSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 924b3ec233SSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status2; 93e26ad7d8SSatish Balay MPI_Comm comm; 944b3ec233SSatish Balay Scalar **rbuf2,**sbuf2; 95e26ad7d8SSatish Balay 96e26ad7d8SSatish Balay PetscFunctionBegin; 97e26ad7d8SSatish Balay comm = C->comm; 98e26ad7d8SSatish Balay tag0 = C->tag; 99e26ad7d8SSatish Balay size = c->size; 100e26ad7d8SSatish Balay rank = c->rank; 101e26ad7d8SSatish Balay m = c->M; 102e26ad7d8SSatish Balay 103e26ad7d8SSatish Balay /* Get some new tags to keep the communication clean */ 104e26ad7d8SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 105e26ad7d8SSatish Balay 106e26ad7d8SSatish Balay /* Check if the col indices are sorted */ 107e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 108888f2ed8SSatish Balay ierr = ISSorted(isrow[i],(PetscTruth*)&j);CHKERRQ(ierr); 109e26ad7d8SSatish Balay if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted"); 110888f2ed8SSatish Balay ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr); 111e26ad7d8SSatish Balay if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted"); 112e26ad7d8SSatish Balay } 113e26ad7d8SSatish Balay 1144b3ec233SSatish Balay len = 2*ismax*(sizeof(int *)+sizeof(int)) + (m+1)*sizeof(int); 115e26ad7d8SSatish Balay irow = (int **)PetscMalloc(len);CHKPTRQ(irow); 116e26ad7d8SSatish Balay icol = irow + ismax; 117e26ad7d8SSatish Balay nrow = (int *) (icol + ismax); 118e26ad7d8SSatish Balay ncol = nrow + ismax; 119e26ad7d8SSatish Balay rtable = ncol + ismax; 120e26ad7d8SSatish Balay 121e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 122e26ad7d8SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 123e26ad7d8SSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 124e26ad7d8SSatish Balay ierr = ISGetSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 125e26ad7d8SSatish Balay ierr = ISGetSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 126e26ad7d8SSatish Balay } 127e26ad7d8SSatish Balay 128e26ad7d8SSatish Balay /* Create hash table for the mapping :row -> proc*/ 129e26ad7d8SSatish Balay for ( i=0,j=0; i<size; i++ ) { 130e26ad7d8SSatish Balay jmax = c->rowners[i+1]; 131e26ad7d8SSatish Balay for ( ; j<jmax; j++ ) { 132e26ad7d8SSatish Balay rtable[j] = i; 133e26ad7d8SSatish Balay } 134e26ad7d8SSatish Balay } 135e26ad7d8SSatish Balay 136e26ad7d8SSatish Balay /* evaluate communication - mesg to who, length of mesg, and buffer space 137e26ad7d8SSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 138e26ad7d8SSatish Balay w1 = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1); /* mesg size */ 139e26ad7d8SSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 140e26ad7d8SSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 141e26ad7d8SSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 142549d3d68SSatish Balay ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/ 143e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 144549d3d68SSatish Balay ierr = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/ 145e26ad7d8SSatish Balay jmax = nrow[i]; 146e26ad7d8SSatish Balay irow_i = irow[i]; 147e26ad7d8SSatish Balay for ( j=0; j<jmax; j++ ) { 148e26ad7d8SSatish Balay row = irow_i[j]; 149e26ad7d8SSatish Balay proc = rtable[row]; 150e26ad7d8SSatish Balay w4[proc]++; 151e26ad7d8SSatish Balay } 152e26ad7d8SSatish Balay for ( j=0; j<size; j++ ) { 153e26ad7d8SSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 154e26ad7d8SSatish Balay } 155e26ad7d8SSatish Balay } 156e26ad7d8SSatish Balay 157e26ad7d8SSatish Balay nrqs = 0; /* no of outgoing messages */ 158e26ad7d8SSatish Balay msz = 0; /* total mesg length (for all procs) */ 159e26ad7d8SSatish Balay w1[rank] = 0; /* no mesg sent to self */ 160e26ad7d8SSatish Balay w3[rank] = 0; 161e26ad7d8SSatish Balay for ( i=0; i<size; i++ ) { 162e26ad7d8SSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 163e26ad7d8SSatish Balay } 164e26ad7d8SSatish Balay pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); /*(proc -array)*/ 165e26ad7d8SSatish Balay for ( i=0, j=0; i<size; i++ ) { 166e26ad7d8SSatish Balay if (w1[i]) { pa[j] = i; j++; } 167e26ad7d8SSatish Balay } 168e26ad7d8SSatish Balay 169e26ad7d8SSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 170e26ad7d8SSatish Balay for ( i=0; i<nrqs; i++ ) { 171e26ad7d8SSatish Balay j = pa[i]; 172e26ad7d8SSatish Balay w1[j] += w2[j] + 2* w3[j]; 173e26ad7d8SSatish Balay msz += w1[j]; 174e26ad7d8SSatish Balay } 175e26ad7d8SSatish Balay /* Do a global reduction to determine how many messages to expect*/ 176e26ad7d8SSatish Balay { 177e26ad7d8SSatish Balay int *rw1, *rw2; 178e26ad7d8SSatish Balay rw1 = (int *)PetscMalloc(2*size*sizeof(int));CHKPTRQ(rw1); 179e26ad7d8SSatish Balay rw2 = rw1+size; 180e26ad7d8SSatish Balay ierr = MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr); 181e26ad7d8SSatish Balay bsz = rw1[rank]; 182e26ad7d8SSatish Balay ierr = MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm);CHKERRQ(ierr); 183e26ad7d8SSatish Balay nrqr = rw2[rank]; 184606d414cSSatish Balay ierr = PetscFree(rw1);CHKERRQ(ierr); 185e26ad7d8SSatish Balay } 186e26ad7d8SSatish Balay 187e26ad7d8SSatish Balay /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 188e26ad7d8SSatish Balay len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int); 189e26ad7d8SSatish Balay rbuf1 = (int**) PetscMalloc(len);CHKPTRQ(rbuf1); 190e26ad7d8SSatish Balay rbuf1[0] = (int *) (rbuf1 + nrqr); 191e26ad7d8SSatish Balay for ( i=1; i<nrqr; ++i ) rbuf1[i] = rbuf1[i-1] + bsz; 192e26ad7d8SSatish Balay 193e26ad7d8SSatish Balay /* Post the receives */ 194e26ad7d8SSatish Balay r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(r_waits1); 195e26ad7d8SSatish Balay for ( i=0; i<nrqr; ++i ) { 196e26ad7d8SSatish Balay ierr = MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr); 197e26ad7d8SSatish Balay } 198e26ad7d8SSatish Balay 199e26ad7d8SSatish Balay /* Allocate Memory for outgoing messages */ 200e26ad7d8SSatish Balay len = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int); 201e26ad7d8SSatish Balay sbuf1 = (int **)PetscMalloc(len);CHKPTRQ(sbuf1); 202e26ad7d8SSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 203549d3d68SSatish Balay ierr = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr); 204e26ad7d8SSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 205e26ad7d8SSatish Balay tmp = (int *) (ptr + size); 206e26ad7d8SSatish Balay ctr = tmp + 2*msz; 207e26ad7d8SSatish Balay 208e26ad7d8SSatish Balay { 209e26ad7d8SSatish Balay int *iptr = tmp,ict = 0; 210e26ad7d8SSatish Balay for ( i=0; i<nrqs; i++ ) { 211e26ad7d8SSatish Balay j = pa[i]; 212e26ad7d8SSatish Balay iptr += ict; 213e26ad7d8SSatish Balay sbuf1[j] = iptr; 214e26ad7d8SSatish Balay ict = w1[j]; 215e26ad7d8SSatish Balay } 216e26ad7d8SSatish Balay } 217e26ad7d8SSatish Balay 218e26ad7d8SSatish Balay /* Form the outgoing messages */ 219e26ad7d8SSatish Balay /* Initialize the header space */ 220e26ad7d8SSatish Balay for ( i=0; i<nrqs; i++ ) { 221e26ad7d8SSatish Balay j = pa[i]; 222e26ad7d8SSatish Balay sbuf1[j][0] = 0; 223549d3d68SSatish Balay ierr = PetscMemzero(sbuf1[j]+1, 2*w3[j]*sizeof(int));CHKERRQ(ierr); 224e26ad7d8SSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 225e26ad7d8SSatish Balay } 226e26ad7d8SSatish Balay 227e26ad7d8SSatish Balay /* Parse the isrow and copy data into outbuf */ 228e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 229549d3d68SSatish Balay ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr); 230e26ad7d8SSatish Balay irow_i = irow[i]; 231e26ad7d8SSatish Balay jmax = nrow[i]; 232e26ad7d8SSatish Balay for ( j=0; j<jmax; j++ ) { /* parse the indices of each IS */ 233e26ad7d8SSatish Balay row = irow_i[j]; 234e26ad7d8SSatish Balay proc = rtable[row]; 235e26ad7d8SSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 236e26ad7d8SSatish Balay ctr[proc]++; 237e26ad7d8SSatish Balay *ptr[proc] = row; 238e26ad7d8SSatish Balay ptr[proc]++; 239e26ad7d8SSatish Balay } 240e26ad7d8SSatish Balay } 241e26ad7d8SSatish Balay /* Update the headers for the current IS */ 242e26ad7d8SSatish Balay for ( j=0; j<size; j++ ) { /* Can Optimise this loop too */ 243e26ad7d8SSatish Balay if ((ctr_j = ctr[j])) { 244e26ad7d8SSatish Balay sbuf1_j = sbuf1[j]; 245e26ad7d8SSatish Balay k = ++sbuf1_j[0]; 246e26ad7d8SSatish Balay sbuf1_j[2*k] = ctr_j; 247e26ad7d8SSatish Balay sbuf1_j[2*k-1] = i; 248e26ad7d8SSatish Balay } 249e26ad7d8SSatish Balay } 250e26ad7d8SSatish Balay } 251e26ad7d8SSatish Balay 252e26ad7d8SSatish Balay /* Now post the sends */ 253e26ad7d8SSatish Balay s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(s_waits1); 254e26ad7d8SSatish Balay for ( i=0; i<nrqs; ++i ) { 255e26ad7d8SSatish Balay j = pa[i]; 256e26ad7d8SSatish Balay ierr = MPI_Isend( sbuf1[j], w1[j], MPI_INT, j, tag0, comm, s_waits1+i);CHKERRQ(ierr); 257e26ad7d8SSatish Balay } 258e26ad7d8SSatish Balay 2594b3ec233SSatish Balay /* Post recieves to capture the row_data from other procs */ 260e26ad7d8SSatish Balay r_waits2 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits2); 2614b3ec233SSatish Balay rbuf2 = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar *));CHKPTRQ(rbuf2); 2624b3ec233SSatish Balay for ( i=0; i<nrqs; i++ ) { 263e26ad7d8SSatish Balay j = pa[i]; 2644b3ec233SSatish Balay count = (w1[j] - (2*sbuf1[j][0] + 1))*N; 2654b3ec233SSatish Balay rbuf2[i] = (Scalar *)PetscMalloc((count+1)*sizeof(Scalar));CHKPTRQ(rbuf2[i]); 2664b3ec233SSatish Balay ierr = MPI_Irecv( rbuf2[i], count, MPIU_SCALAR, j, tag1, comm, r_waits2+i);CHKERRQ(ierr); 267e26ad7d8SSatish Balay } 268e26ad7d8SSatish Balay 2694b3ec233SSatish Balay /* Receive messages(row_nos) and then, pack and send off the rowvalues 2704b3ec233SSatish Balay to the correct processors */ 271e26ad7d8SSatish Balay 272e26ad7d8SSatish Balay s_waits2 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits2); 273e26ad7d8SSatish Balay r_status1 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(r_status1); 2744b3ec233SSatish Balay sbuf2 = (Scalar **) PetscMalloc((nrqr+1)*sizeof(Scalar*));CHKPTRQ(sbuf2); 275e26ad7d8SSatish Balay 276e26ad7d8SSatish Balay { 2774b3ec233SSatish Balay Scalar *sbuf2_i,*v_start; 2784b3ec233SSatish Balay int s_proc; 279e26ad7d8SSatish Balay for ( i=0; i<nrqr; ++i ) { 280e26ad7d8SSatish Balay ierr = MPI_Waitany(nrqr, r_waits1, &index, r_status1+i);CHKERRQ(ierr); 2814b3ec233SSatish Balay s_proc = r_status1[i].MPI_SOURCE; /* send processor */ 2824b3ec233SSatish Balay rbuf1_i = rbuf1[index]; /* Actual message from s_proc */ 2834b3ec233SSatish Balay /* no of rows = end - start; since start is array index[], 0index, whel end 2844b3ec233SSatish Balay is length of the buffer - which is 1index */ 285e26ad7d8SSatish Balay start = 2*rbuf1_i[0] + 1; 286e26ad7d8SSatish Balay MPI_Get_count(r_status1+i,MPI_INT, &end); 2874b3ec233SSatish Balay /* allocate memory sufficinet to hold all the row values */ 2884b3ec233SSatish Balay sbuf2[index] = (Scalar *)PetscMalloc((end-start)*N*sizeof(Scalar));CHKPTRQ(sbuf2[index]); 289e26ad7d8SSatish Balay sbuf2_i = sbuf2[index]; 2904b3ec233SSatish Balay /* Now pack the data */ 291e26ad7d8SSatish Balay for ( j=start; j<end; j++ ) { 2924b3ec233SSatish Balay row = rbuf1_i[j] - rstart; 2934b3ec233SSatish Balay v_start = a->v + row; 2944b3ec233SSatish Balay for ( k=0; k<N; k++ ) { 2954b3ec233SSatish Balay sbuf2_i[0] = v_start[0]; 2964b3ec233SSatish Balay sbuf2_i++; v_start+=a->m; 297e26ad7d8SSatish Balay } 298e26ad7d8SSatish Balay } 2994b3ec233SSatish Balay /* Now send off the data */ 3004b3ec233SSatish Balay ierr = MPI_Isend(sbuf2[index],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);CHKERRQ(ierr); 3014b3ec233SSatish Balay } 3024b3ec233SSatish Balay } 3034b3ec233SSatish Balay /* End Send-Recv of IS + row_numbers */ 304606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 305606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 306e26ad7d8SSatish Balay s_status1 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(s_status1); 307e26ad7d8SSatish Balay ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr); 308606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 309606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 310e26ad7d8SSatish Balay 3114b3ec233SSatish Balay /* Create the submatrices */ 3124b3ec233SSatish Balay if (scall == MAT_REUSE_MATRIX) { 313e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 3144b3ec233SSatish Balay mat = (Mat_SeqDense *)(submats[i]->data); 3154b3ec233SSatish Balay if ((mat->m != nrow[i]) || (mat->n != ncol[i])) { 3164b3ec233SSatish Balay SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size"); 317e26ad7d8SSatish Balay } 318549d3d68SSatish Balay ierr = PetscMemzero(mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr); 3194b3ec233SSatish Balay submats[i]->factor = C->factor; 320e26ad7d8SSatish Balay } 3214b3ec233SSatish Balay } else { 322e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 3234b3ec233SSatish Balay ierr = MatCreateSeqDense(PETSC_COMM_SELF,nrow[i],ncol[i],PETSC_NULL,submats+i);CHKERRQ(ierr); 3244b3ec233SSatish Balay } 3254b3ec233SSatish Balay } 3264b3ec233SSatish Balay 3274b3ec233SSatish Balay /* Assemble the matrices */ 3284b3ec233SSatish Balay { 3294b3ec233SSatish Balay int col; 3304b3ec233SSatish Balay Scalar *imat_v,*mat_v,*imat_vi,*mat_vi; 3314b3ec233SSatish Balay 3324b3ec233SSatish Balay for ( i=0; i<ismax; i++ ) { 3334b3ec233SSatish Balay mat = (Mat_SeqDense *) submats[i]->data; 3344b3ec233SSatish Balay mat_v = a->v; 3354b3ec233SSatish Balay imat_v = mat->v; 336e26ad7d8SSatish Balay irow_i = irow[i]; 3374b3ec233SSatish Balay m = nrow[i]; 3384b3ec233SSatish Balay for ( j=0; j<m; j++ ) { 339e26ad7d8SSatish Balay row = irow_i[j] ; 340e26ad7d8SSatish Balay proc = rtable[row]; 341e26ad7d8SSatish Balay if (proc == rank) { 3424b3ec233SSatish Balay row = row - rstart; 3434b3ec233SSatish Balay mat_vi = mat_v + row; 3444b3ec233SSatish Balay imat_vi = imat_v + j; 3454b3ec233SSatish Balay for ( k=0; k<ncol[i]; k++ ) { 3464b3ec233SSatish Balay col = icol[i][k]; 3474b3ec233SSatish Balay imat_vi[k*m] = mat_vi[col*a->m]; 348e26ad7d8SSatish Balay } 3494b3ec233SSatish Balay } 350e26ad7d8SSatish Balay } 351e26ad7d8SSatish Balay } 352e26ad7d8SSatish Balay } 353e26ad7d8SSatish Balay 3544b3ec233SSatish Balay /* Create row map. This maps c->row to submat->row for each submat*/ 3554b3ec233SSatish Balay /* this is a very expensive operation wrt memory usage */ 356e26ad7d8SSatish Balay len = (1+ismax)*sizeof(int *) + ismax*c->M*sizeof(int); 357e26ad7d8SSatish Balay rmap = (int **)PetscMalloc(len);CHKPTRQ(rmap); 358e26ad7d8SSatish Balay rmap[0] = (int *)(rmap + ismax); 359549d3d68SSatish Balay ierr = PetscMemzero(rmap[0],ismax*c->M*sizeof(int));CHKERRQ(ierr); 360e26ad7d8SSatish Balay for ( i=1; i<ismax; i++ ) { rmap[i] = rmap[i-1] + c->M;} 361e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 362e26ad7d8SSatish Balay rmap_i = rmap[i]; 363e26ad7d8SSatish Balay irow_i = irow[i]; 364e26ad7d8SSatish Balay jmax = nrow[i]; 365e26ad7d8SSatish Balay for ( j=0; j<jmax; j++ ) { 366e26ad7d8SSatish Balay rmap_i[irow_i[j]] = j; 367e26ad7d8SSatish Balay } 368e26ad7d8SSatish Balay } 369e26ad7d8SSatish Balay 3704b3ec233SSatish Balay /* Now Receive the row_values and assemble the rest of the matrix */ 371e26ad7d8SSatish Balay 3724b3ec233SSatish Balay r_status2 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status2); 3734b3ec233SSatish Balay 3744b3ec233SSatish Balay { 3754b3ec233SSatish Balay int is_max,tmp1,col,*sbuf1_i,is_sz; 3764b3ec233SSatish Balay Scalar *rbuf2_i,*imat_v,*imat_vi; 3774b3ec233SSatish Balay 3784b3ec233SSatish Balay for ( tmp1=0; tmp1<nrqs; tmp1++ ) { /* For each message */ 3794b3ec233SSatish Balay ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr); 3804b3ec233SSatish Balay /* Now dig out the corresponding sbuf1, which contains the IS data_structure */ 3814b3ec233SSatish Balay sbuf1_i = sbuf1[pa[i]]; 3824b3ec233SSatish Balay is_max = sbuf1_i[0]; 3834b3ec233SSatish Balay ct1 = 2*is_max+1; 384e26ad7d8SSatish Balay rbuf2_i = rbuf2[i]; 3854b3ec233SSatish Balay for ( j=1; j<=is_max; j++ ) { /* For each IS belonging to the message */ 386e26ad7d8SSatish Balay is_no = sbuf1_i[2*j-1]; 3874b3ec233SSatish Balay is_sz = sbuf1_i[2*j]; 388e26ad7d8SSatish Balay mat = (Mat_SeqDense *) submats[is_no]->data; 3894b3ec233SSatish Balay imat_v = mat->v; 3904b3ec233SSatish Balay rmap_i = rmap[is_no]; 3914b3ec233SSatish Balay m = nrow[is_no]; 3924b3ec233SSatish Balay for ( k=0; k<is_sz; k++,rbuf2_i+=N) { /* For each row */ 3934b3ec233SSatish Balay row = sbuf1_i[ct1]; ct1++; 394e26ad7d8SSatish Balay row = rmap_i[row]; 3954b3ec233SSatish Balay imat_vi = imat_v + row; 3964b3ec233SSatish Balay for ( l=0; l<ncol[is_no]; l++ ) { /* For each col */ 3974b3ec233SSatish Balay col = icol[is_no][l]; 3984b3ec233SSatish Balay imat_vi[l*m] = rbuf2_i[col]; 399e26ad7d8SSatish Balay } 400e26ad7d8SSatish Balay } 401e26ad7d8SSatish Balay } 402e26ad7d8SSatish Balay } 4034b3ec233SSatish Balay } 4044b3ec233SSatish Balay /* End Send-Recv of row_values */ 405606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 406606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 4074b3ec233SSatish Balay s_status2 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status1); 4084b3ec233SSatish Balay ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr); 409606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 410606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 411e26ad7d8SSatish Balay 412e26ad7d8SSatish Balay /* Restore the indices */ 413e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 414e26ad7d8SSatish Balay ierr = ISRestoreIndices(isrow[i], irow+i);CHKERRQ(ierr); 415e26ad7d8SSatish Balay ierr = ISRestoreIndices(iscol[i], icol+i);CHKERRQ(ierr); 416e26ad7d8SSatish Balay } 417e26ad7d8SSatish Balay 418e26ad7d8SSatish Balay /* Destroy allocated memory */ 419606d414cSSatish Balay ierr = PetscFree(irow);CHKERRQ(ierr); 420606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 421606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 422e26ad7d8SSatish Balay 4234b3ec233SSatish Balay 4244b3ec233SSatish Balay for ( i=0; i<nrqs; ++i ) { 425606d414cSSatish Balay ierr = PetscFree(rbuf2[i]);CHKERRQ(ierr); 4264b3ec233SSatish Balay } 427606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 428606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 429606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 4304b3ec233SSatish Balay 431e26ad7d8SSatish Balay for ( i=0; i<nrqr; ++i ) { 432606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 433e26ad7d8SSatish Balay } 434e26ad7d8SSatish Balay 435606d414cSSatish Balay ierr = PetscFree(sbuf2);CHKERRQ(ierr); 436606d414cSSatish Balay ierr = PetscFree(rmap);CHKERRQ(ierr); 437e26ad7d8SSatish Balay 438e26ad7d8SSatish Balay ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 439e26ad7d8SSatish Balay 440e26ad7d8SSatish Balay for ( i=0; i<ismax; i++ ) { 441e26ad7d8SSatish Balay ierr = MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 442e26ad7d8SSatish Balay ierr = MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 443e26ad7d8SSatish Balay } 444e26ad7d8SSatish Balay 445e26ad7d8SSatish Balay PetscFunctionReturn(0); 446e26ad7d8SSatish Balay } 447e26ad7d8SSatish Balay 448