xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision 606d414c6e034e02e67059b83ebaefc3ebe99698)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*606d414cSSatish Balay static char vcid[] = "$Id: mmdense.c,v 1.21 1999/06/08 22:55:41 balay Exp balay $";
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 
30426b7fd2SBarry Smith   /* the partitioning of this vector is not right! you don't really know
31426b7fd2SBarry Smith      what it is; this will only work for default (PETSC_DECIDE) matrix layouts*/
32426b7fd2SBarry Smith   ierr = VecCreateMPI(mat->comm,PETSC_DECIDE,mdn->N,&gvec);CHKERRQ(ierr);
33e19c7942SLois Curfman McInnes 
340500aeceSBarry Smith   /* Generate the scatter context */
35d9f96c7cSLois Curfman McInnes   ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx);CHKERRQ(ierr);
360500aeceSBarry Smith   PLogObjectParent(mat,mdn->Mvctx);
370500aeceSBarry Smith   PLogObjectParent(mat,mdn->lvec);
38d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,tofrom);
39d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,gvec);
40e19c7942SLois Curfman McInnes 
41d9f96c7cSLois Curfman McInnes   ierr = ISDestroy(tofrom);CHKERRQ(ierr);
42d9f96c7cSLois Curfman McInnes   ierr = VecDestroy(gvec);CHKERRQ(ierr);
433a40ed3dSBarry Smith   PetscFunctionReturn(0);
440500aeceSBarry Smith }
450500aeceSBarry Smith 
467b2a1423SBarry Smith extern int MatGetSubMatrices_MPIDense_Local(Mat,int,IS*,IS*,MatReuse,Mat*);
47e26ad7d8SSatish Balay #undef __FUNC__
48e26ad7d8SSatish Balay #define __FUNC__ "MatGetSubMatrices_MPIDense"
49e26ad7d8SSatish Balay int MatGetSubMatrices_MPIDense(Mat C,int ismax,IS *isrow,IS *iscol,
507b2a1423SBarry Smith                              MatReuse scall,Mat **submat)
51e26ad7d8SSatish Balay {
52e26ad7d8SSatish Balay   Mat_MPIDense  *c = (Mat_MPIDense *) C->data;
53e26ad7d8SSatish Balay   int           nmax,nstages_local,nstages,i,pos,max_no,ierr;
54e26ad7d8SSatish Balay 
55e26ad7d8SSatish Balay   PetscFunctionBegin;
56e26ad7d8SSatish Balay   /* Allocate memory to hold all the submatrices */
57e26ad7d8SSatish Balay   if (scall != MAT_REUSE_MATRIX) {
58e26ad7d8SSatish Balay     *submat = (Mat *)PetscMalloc((ismax+1)*sizeof(Mat));CHKPTRQ(*submat);
59e26ad7d8SSatish Balay   }
60e26ad7d8SSatish Balay   /* Determine the number of stages through which submatrices are done */
61e26ad7d8SSatish Balay   nmax          = 20*1000000 / (c->N * sizeof(int));
62e26ad7d8SSatish Balay   if (nmax == 0) nmax = 1;
63e26ad7d8SSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
64e26ad7d8SSatish Balay 
65e26ad7d8SSatish Balay   /* Make sure every processor loops through the nstages */
66e26ad7d8SSatish Balay   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr);
670500aeceSBarry Smith 
680500aeceSBarry Smith 
69e26ad7d8SSatish Balay   for ( i=0,pos=0; i<nstages; i++ ) {
70e26ad7d8SSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
71e26ad7d8SSatish Balay     else if (pos == ismax) max_no = 0;
72e26ad7d8SSatish Balay     else                   max_no = ismax-pos;
73e26ad7d8SSatish Balay     ierr = MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);CHKERRQ(ierr);
74e26ad7d8SSatish Balay     pos += max_no;
75e26ad7d8SSatish Balay   }
76e26ad7d8SSatish Balay   PetscFunctionReturn(0);
77e26ad7d8SSatish Balay }
78e26ad7d8SSatish Balay /* -------------------------------------------------------------------------*/
79e26ad7d8SSatish Balay #undef __FUNC__
80e26ad7d8SSatish Balay #define __FUNC__ "MatGetSubMatrices_MPIDense_Local"
817b2a1423SBarry Smith int MatGetSubMatrices_MPIDense_Local(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat *submats)
82e26ad7d8SSatish Balay {
83e26ad7d8SSatish Balay   Mat_MPIDense  *c = (Mat_MPIDense *) C->data;
84e26ad7d8SSatish Balay   Mat            A = c->A;
85e26ad7d8SSatish Balay   Mat_SeqDense  *a = (Mat_SeqDense*)A->data, *mat;
864b3ec233SSatish Balay   int           N = c->N, rstart = c->rstart,count;
87e26ad7d8SSatish Balay   int           **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size;
884b3ec233SSatish Balay   int           **sbuf1, rank, m,i,j,k,l,ct1,ierr, **rbuf1,row,proc;
894b3ec233SSatish Balay   int           nrqs, msz, **ptr,index,*ctr,*pa,*tmp,bsz,nrqr;
904b3ec233SSatish Balay   int           is_no,jmax,*irow_i,**rmap,*rmap_i;
914b3ec233SSatish Balay   int           len,ctr_j,*sbuf1_j,*rbuf1_i;
924b3ec233SSatish Balay   int           tag0,tag1;
934b3ec233SSatish Balay   MPI_Request   *s_waits1,*r_waits1,*s_waits2,*r_waits2;
944b3ec233SSatish Balay   MPI_Status    *r_status1,*r_status2,*s_status1,*s_status2;
95e26ad7d8SSatish Balay   MPI_Comm      comm;
964b3ec233SSatish Balay   Scalar        **rbuf2,**sbuf2;
97e26ad7d8SSatish Balay 
98e26ad7d8SSatish Balay   PetscFunctionBegin;
99e26ad7d8SSatish Balay   comm   = C->comm;
100e26ad7d8SSatish Balay   tag0   = C->tag;
101e26ad7d8SSatish Balay   size   = c->size;
102e26ad7d8SSatish Balay   rank   = c->rank;
103e26ad7d8SSatish Balay   m      = c->M;
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++ ) {
110888f2ed8SSatish Balay     ierr = ISSorted(isrow[i],(PetscTruth*)&j);CHKERRQ(ierr);
111e26ad7d8SSatish Balay     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted");
112888f2ed8SSatish Balay     ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr);
113e26ad7d8SSatish Balay     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted");
114e26ad7d8SSatish Balay   }
115e26ad7d8SSatish Balay 
1164b3ec233SSatish Balay   len    =  2*ismax*(sizeof(int *)+sizeof(int)) + (m+1)*sizeof(int);
117e26ad7d8SSatish Balay   irow   = (int **)PetscMalloc(len);CHKPTRQ(irow);
118e26ad7d8SSatish Balay   icol   = irow + ismax;
119e26ad7d8SSatish Balay   nrow   = (int *) (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);
126e26ad7d8SSatish Balay     ierr = ISGetSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
127e26ad7d8SSatish Balay     ierr = ISGetSize(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++ ) {
132e26ad7d8SSatish Balay     jmax = c->rowners[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*/
140e26ad7d8SSatish Balay   w1     = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1); /* mesg size */
141e26ad7d8SSatish Balay   w2     = w1 + size;      /* if w2[i] marked, then a message to proc i*/
142e26ad7d8SSatish Balay   w3     = w2 + size;      /* no of IS that needs to be sent to proc i */
143e26ad7d8SSatish Balay   w4     = w3 + size;      /* temp work space used in determining w1, w2, w3 */
144549d3d68SSatish Balay   ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/
145e26ad7d8SSatish Balay   for ( i=0; i<ismax; i++ ) {
146549d3d68SSatish Balay     ierr   = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/
147e26ad7d8SSatish Balay     jmax   = nrow[i];
148e26ad7d8SSatish Balay     irow_i = irow[i];
149e26ad7d8SSatish Balay     for ( j=0; j<jmax; j++ ) {
150e26ad7d8SSatish Balay       row  = irow_i[j];
151e26ad7d8SSatish Balay       proc = rtable[row];
152e26ad7d8SSatish Balay       w4[proc]++;
153e26ad7d8SSatish Balay     }
154e26ad7d8SSatish Balay     for ( j=0; j<size; j++ ) {
155e26ad7d8SSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
156e26ad7d8SSatish Balay     }
157e26ad7d8SSatish Balay   }
158e26ad7d8SSatish Balay 
159e26ad7d8SSatish Balay   nrqs     = 0;              /* no of outgoing messages */
160e26ad7d8SSatish Balay   msz      = 0;              /* total mesg length (for all procs) */
161e26ad7d8SSatish Balay   w1[rank] = 0;              /* no mesg sent to self */
162e26ad7d8SSatish Balay   w3[rank] = 0;
163e26ad7d8SSatish Balay   for ( i=0; i<size; i++ ) {
164e26ad7d8SSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
165e26ad7d8SSatish Balay   }
166e26ad7d8SSatish Balay   pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); /*(proc -array)*/
167e26ad7d8SSatish Balay   for ( i=0, j=0; i<size; i++ ) {
168e26ad7d8SSatish Balay     if (w1[i]) { pa[j] = i; j++; }
169e26ad7d8SSatish Balay   }
170e26ad7d8SSatish Balay 
171e26ad7d8SSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
172e26ad7d8SSatish Balay   for ( i=0; i<nrqs; i++ ) {
173e26ad7d8SSatish Balay     j     = pa[i];
174e26ad7d8SSatish Balay     w1[j] += w2[j] + 2* w3[j];
175e26ad7d8SSatish Balay     msz   += w1[j];
176e26ad7d8SSatish Balay   }
177e26ad7d8SSatish Balay   /* Do a global reduction to determine how many messages to expect*/
178e26ad7d8SSatish Balay   {
179e26ad7d8SSatish Balay     int *rw1, *rw2;
180e26ad7d8SSatish Balay     rw1   = (int *)PetscMalloc(2*size*sizeof(int));CHKPTRQ(rw1);
181e26ad7d8SSatish Balay     rw2   = rw1+size;
182e26ad7d8SSatish Balay     ierr  = MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
183e26ad7d8SSatish Balay     bsz   = rw1[rank];
184e26ad7d8SSatish Balay     ierr  = MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm);CHKERRQ(ierr);
185e26ad7d8SSatish Balay     nrqr  = rw2[rank];
186*606d414cSSatish Balay     ierr  = PetscFree(rw1);CHKERRQ(ierr);
187e26ad7d8SSatish Balay   }
188e26ad7d8SSatish Balay 
189e26ad7d8SSatish Balay   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
190e26ad7d8SSatish Balay   len      = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
191e26ad7d8SSatish Balay   rbuf1    = (int**) PetscMalloc(len);CHKPTRQ(rbuf1);
192e26ad7d8SSatish Balay   rbuf1[0] = (int *) (rbuf1 + nrqr);
193e26ad7d8SSatish Balay   for ( i=1; i<nrqr; ++i ) rbuf1[i] = rbuf1[i-1] + bsz;
194e26ad7d8SSatish Balay 
195e26ad7d8SSatish Balay   /* Post the receives */
196e26ad7d8SSatish Balay   r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(r_waits1);
197e26ad7d8SSatish Balay   for ( i=0; i<nrqr; ++i ) {
198e26ad7d8SSatish Balay     ierr = MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr);
199e26ad7d8SSatish Balay   }
200e26ad7d8SSatish Balay 
201e26ad7d8SSatish Balay   /* Allocate Memory for outgoing messages */
202e26ad7d8SSatish Balay   len   = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int);
203e26ad7d8SSatish Balay   sbuf1 = (int **)PetscMalloc(len);CHKPTRQ(sbuf1);
204e26ad7d8SSatish Balay   ptr   = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
205549d3d68SSatish Balay   ierr  = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr);
206e26ad7d8SSatish Balay   /* allocate memory for outgoing data + buf to receive the first reply */
207e26ad7d8SSatish Balay   tmp   = (int *) (ptr + size);
208e26ad7d8SSatish Balay   ctr   = tmp + 2*msz;
209e26ad7d8SSatish Balay 
210e26ad7d8SSatish Balay   {
211e26ad7d8SSatish Balay     int *iptr = tmp,ict = 0;
212e26ad7d8SSatish Balay     for ( i=0; i<nrqs; i++ ) {
213e26ad7d8SSatish Balay       j         = pa[i];
214e26ad7d8SSatish Balay       iptr     += ict;
215e26ad7d8SSatish Balay       sbuf1[j]  = iptr;
216e26ad7d8SSatish Balay       ict       = w1[j];
217e26ad7d8SSatish Balay     }
218e26ad7d8SSatish Balay   }
219e26ad7d8SSatish Balay 
220e26ad7d8SSatish Balay   /* Form the outgoing messages */
221e26ad7d8SSatish Balay   /* Initialize the header space */
222e26ad7d8SSatish Balay   for ( i=0; i<nrqs; i++ ) {
223e26ad7d8SSatish Balay     j           = pa[i];
224e26ad7d8SSatish Balay     sbuf1[j][0] = 0;
225549d3d68SSatish Balay     ierr        = PetscMemzero(sbuf1[j]+1, 2*w3[j]*sizeof(int));CHKERRQ(ierr);
226e26ad7d8SSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
227e26ad7d8SSatish Balay   }
228e26ad7d8SSatish Balay 
229e26ad7d8SSatish Balay   /* Parse the isrow and copy data into outbuf */
230e26ad7d8SSatish Balay   for ( i=0; i<ismax; i++ ) {
231549d3d68SSatish Balay     ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr);
232e26ad7d8SSatish Balay     irow_i = irow[i];
233e26ad7d8SSatish Balay     jmax   = nrow[i];
234e26ad7d8SSatish Balay     for ( j=0; j<jmax; j++ ) {  /* parse the indices of each IS */
235e26ad7d8SSatish Balay       row  = irow_i[j];
236e26ad7d8SSatish Balay       proc = rtable[row];
237e26ad7d8SSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
238e26ad7d8SSatish Balay         ctr[proc]++;
239e26ad7d8SSatish Balay         *ptr[proc] = row;
240e26ad7d8SSatish Balay         ptr[proc]++;
241e26ad7d8SSatish Balay       }
242e26ad7d8SSatish Balay     }
243e26ad7d8SSatish Balay     /* Update the headers for the current IS */
244e26ad7d8SSatish Balay     for ( j=0; j<size; j++ ) { /* Can Optimise this loop too */
245e26ad7d8SSatish Balay       if ((ctr_j = ctr[j])) {
246e26ad7d8SSatish Balay         sbuf1_j        = sbuf1[j];
247e26ad7d8SSatish Balay         k              = ++sbuf1_j[0];
248e26ad7d8SSatish Balay         sbuf1_j[2*k]   = ctr_j;
249e26ad7d8SSatish Balay         sbuf1_j[2*k-1] = i;
250e26ad7d8SSatish Balay       }
251e26ad7d8SSatish Balay     }
252e26ad7d8SSatish Balay   }
253e26ad7d8SSatish Balay 
254e26ad7d8SSatish Balay   /*  Now  post the sends */
255e26ad7d8SSatish Balay   s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(s_waits1);
256e26ad7d8SSatish Balay   for ( i=0; i<nrqs; ++i ) {
257e26ad7d8SSatish Balay     j = pa[i];
258e26ad7d8SSatish Balay     ierr = MPI_Isend( sbuf1[j], w1[j], MPI_INT, j, tag0, comm, s_waits1+i);CHKERRQ(ierr);
259e26ad7d8SSatish Balay   }
260e26ad7d8SSatish Balay 
2614b3ec233SSatish Balay   /* Post recieves to capture the row_data from other procs */
262e26ad7d8SSatish Balay   r_waits2 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits2);
2634b3ec233SSatish Balay   rbuf2    = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar *));CHKPTRQ(rbuf2);
2644b3ec233SSatish Balay   for ( i=0; i<nrqs; i++ ) {
265e26ad7d8SSatish Balay     j        = pa[i];
2664b3ec233SSatish Balay     count    = (w1[j] - (2*sbuf1[j][0] + 1))*N;
2674b3ec233SSatish Balay     rbuf2[i] = (Scalar *)PetscMalloc((count+1)*sizeof(Scalar));CHKPTRQ(rbuf2[i]);
2684b3ec233SSatish Balay     ierr     = MPI_Irecv( rbuf2[i], count, MPIU_SCALAR, j, tag1, comm, r_waits2+i);CHKERRQ(ierr);
269e26ad7d8SSatish Balay   }
270e26ad7d8SSatish Balay 
2714b3ec233SSatish Balay   /* Receive messages(row_nos) and then, pack and send off the rowvalues
2724b3ec233SSatish Balay      to the correct processors */
273e26ad7d8SSatish Balay 
274e26ad7d8SSatish Balay   s_waits2  = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits2);
275e26ad7d8SSatish Balay   r_status1 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(r_status1);
2764b3ec233SSatish Balay   sbuf2       = (Scalar **) PetscMalloc((nrqr+1)*sizeof(Scalar*));CHKPTRQ(sbuf2);
277e26ad7d8SSatish Balay 
278e26ad7d8SSatish Balay   {
2794b3ec233SSatish Balay     Scalar *sbuf2_i,*v_start;
2804b3ec233SSatish Balay     int    s_proc;
281e26ad7d8SSatish Balay     for ( i=0; i<nrqr; ++i ) {
282e26ad7d8SSatish Balay       ierr = MPI_Waitany(nrqr, r_waits1, &index, r_status1+i);CHKERRQ(ierr);
2834b3ec233SSatish Balay       s_proc          = r_status1[i].MPI_SOURCE; /* send processor */
2844b3ec233SSatish Balay       rbuf1_i         = rbuf1[index]; /* Actual message from s_proc */
2854b3ec233SSatish Balay       /* no of rows = end - start; since start is array index[], 0index, whel end
2864b3ec233SSatish Balay          is length of the buffer - which is 1index */
287e26ad7d8SSatish Balay       start           = 2*rbuf1_i[0] + 1;
288e26ad7d8SSatish Balay       MPI_Get_count(r_status1+i,MPI_INT, &end);
2894b3ec233SSatish Balay       /* allocate memory sufficinet to hold all the row values */
2904b3ec233SSatish Balay       sbuf2[index] = (Scalar *)PetscMalloc((end-start)*N*sizeof(Scalar));CHKPTRQ(sbuf2[index]);
291e26ad7d8SSatish Balay       sbuf2_i      = sbuf2[index];
2924b3ec233SSatish Balay       /* Now pack the data */
293e26ad7d8SSatish Balay       for ( j=start; j<end; j++ ) {
2944b3ec233SSatish Balay         row = rbuf1_i[j] - rstart;
2954b3ec233SSatish Balay         v_start = a->v + row;
2964b3ec233SSatish Balay         for ( k=0; k<N; k++ ) {
2974b3ec233SSatish Balay           sbuf2_i[0] = v_start[0];
2984b3ec233SSatish Balay           sbuf2_i++; v_start+=a->m;
299e26ad7d8SSatish Balay         }
300e26ad7d8SSatish Balay       }
3014b3ec233SSatish Balay       /* Now send off the data */
3024b3ec233SSatish Balay       ierr = MPI_Isend(sbuf2[index],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);CHKERRQ(ierr);
3034b3ec233SSatish Balay     }
3044b3ec233SSatish Balay   }
3054b3ec233SSatish Balay   /* End Send-Recv of IS + row_numbers */
306*606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
307*606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
308e26ad7d8SSatish Balay   s_status1 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(s_status1);
309e26ad7d8SSatish Balay   ierr      = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);
310*606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
311*606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
312e26ad7d8SSatish Balay 
3134b3ec233SSatish Balay   /* Create the submatrices */
3144b3ec233SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
315e26ad7d8SSatish Balay     for ( i=0; i<ismax; i++ ) {
3164b3ec233SSatish Balay       mat = (Mat_SeqDense *)(submats[i]->data);
3174b3ec233SSatish Balay       if ((mat->m != nrow[i]) || (mat->n != ncol[i])) {
3184b3ec233SSatish Balay         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
319e26ad7d8SSatish Balay       }
320549d3d68SSatish Balay       ierr = PetscMemzero(mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr);
3214b3ec233SSatish Balay       submats[i]->factor = C->factor;
322e26ad7d8SSatish Balay     }
3234b3ec233SSatish Balay   } else {
324e26ad7d8SSatish Balay     for ( i=0; i<ismax; i++ ) {
3254b3ec233SSatish Balay       ierr = MatCreateSeqDense(PETSC_COMM_SELF,nrow[i],ncol[i],PETSC_NULL,submats+i);CHKERRQ(ierr);
3264b3ec233SSatish Balay     }
3274b3ec233SSatish Balay   }
3284b3ec233SSatish Balay 
3294b3ec233SSatish Balay   /* Assemble the matrices */
3304b3ec233SSatish Balay   {
3314b3ec233SSatish Balay     int    col;
3324b3ec233SSatish Balay     Scalar *imat_v,*mat_v,*imat_vi,*mat_vi;
3334b3ec233SSatish Balay 
3344b3ec233SSatish Balay     for ( i=0; i<ismax; i++ ) {
3354b3ec233SSatish Balay       mat       = (Mat_SeqDense *) submats[i]->data;
3364b3ec233SSatish Balay       mat_v     = a->v;
3374b3ec233SSatish Balay       imat_v    = mat->v;
338e26ad7d8SSatish Balay       irow_i    = irow[i];
3394b3ec233SSatish Balay       m         = nrow[i];
3404b3ec233SSatish Balay       for ( j=0; j<m; j++ ) {
341e26ad7d8SSatish Balay         row      = irow_i[j] ;
342e26ad7d8SSatish Balay         proc     = rtable[row];
343e26ad7d8SSatish Balay         if (proc == rank) {
3444b3ec233SSatish Balay           row      = row - rstart;
3454b3ec233SSatish Balay           mat_vi   = mat_v + row;
3464b3ec233SSatish Balay           imat_vi  = imat_v + j;
3474b3ec233SSatish Balay           for ( k=0; k<ncol[i]; k++ ) {
3484b3ec233SSatish Balay             col = icol[i][k];
3494b3ec233SSatish Balay             imat_vi[k*m] = mat_vi[col*a->m];
350e26ad7d8SSatish Balay           }
3514b3ec233SSatish Balay         }
352e26ad7d8SSatish Balay       }
353e26ad7d8SSatish Balay     }
354e26ad7d8SSatish Balay   }
355e26ad7d8SSatish Balay 
3564b3ec233SSatish Balay   /* Create row map. This maps c->row to submat->row for each submat*/
3574b3ec233SSatish Balay   /* this is a very expensive operation wrt memory usage */
358e26ad7d8SSatish Balay   len     = (1+ismax)*sizeof(int *) + ismax*c->M*sizeof(int);
359e26ad7d8SSatish Balay   rmap    = (int **)PetscMalloc(len);CHKPTRQ(rmap);
360e26ad7d8SSatish Balay   rmap[0] = (int *)(rmap + ismax);
361549d3d68SSatish Balay   ierr    = PetscMemzero(rmap[0],ismax*c->M*sizeof(int));CHKERRQ(ierr);
362e26ad7d8SSatish Balay   for ( i=1; i<ismax; i++ ) { rmap[i] = rmap[i-1] + c->M;}
363e26ad7d8SSatish Balay   for ( i=0; i<ismax; i++ ) {
364e26ad7d8SSatish Balay     rmap_i = rmap[i];
365e26ad7d8SSatish Balay     irow_i = irow[i];
366e26ad7d8SSatish Balay     jmax   = nrow[i];
367e26ad7d8SSatish Balay     for ( j=0; j<jmax; j++ ) {
368e26ad7d8SSatish Balay       rmap_i[irow_i[j]] = j;
369e26ad7d8SSatish Balay     }
370e26ad7d8SSatish Balay   }
371e26ad7d8SSatish Balay 
3724b3ec233SSatish Balay   /* Now Receive the row_values and assemble the rest of the matrix */
373e26ad7d8SSatish Balay 
3744b3ec233SSatish Balay   r_status2 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status2);
3754b3ec233SSatish Balay 
3764b3ec233SSatish Balay   {
3774b3ec233SSatish Balay     int    is_max,tmp1,col,*sbuf1_i,is_sz;
3784b3ec233SSatish Balay     Scalar *rbuf2_i,*imat_v,*imat_vi;
3794b3ec233SSatish Balay 
3804b3ec233SSatish Balay     for ( tmp1=0; tmp1<nrqs; tmp1++ ) { /* For each message */
3814b3ec233SSatish Balay       ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr);
3824b3ec233SSatish Balay       /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
3834b3ec233SSatish Balay       sbuf1_i = sbuf1[pa[i]];
3844b3ec233SSatish Balay       is_max  = sbuf1_i[0];
3854b3ec233SSatish Balay       ct1     = 2*is_max+1;
386e26ad7d8SSatish Balay       rbuf2_i = rbuf2[i];
3874b3ec233SSatish Balay       for ( j=1; j<=is_max; j++ ) { /* For each IS belonging to the message */
388e26ad7d8SSatish Balay         is_no     = sbuf1_i[2*j-1];
3894b3ec233SSatish Balay         is_sz     = sbuf1_i[2*j];
390e26ad7d8SSatish Balay         mat       = (Mat_SeqDense *) submats[is_no]->data;
3914b3ec233SSatish Balay         imat_v    = mat->v;
3924b3ec233SSatish Balay         rmap_i    = rmap[is_no];
3934b3ec233SSatish Balay         m         = nrow[is_no];
3944b3ec233SSatish Balay         for ( k=0; k<is_sz; k++,rbuf2_i+=N) {  /* For each row */
3954b3ec233SSatish Balay           row      = sbuf1_i[ct1]; ct1++;
396e26ad7d8SSatish Balay           row      = rmap_i[row];
3974b3ec233SSatish Balay           imat_vi  = imat_v + row;
3984b3ec233SSatish Balay           for ( l=0; l<ncol[is_no]; l++ ) { /* For each col */
3994b3ec233SSatish Balay             col = icol[is_no][l];
4004b3ec233SSatish Balay             imat_vi[l*m] = rbuf2_i[col];
401e26ad7d8SSatish Balay           }
402e26ad7d8SSatish Balay         }
403e26ad7d8SSatish Balay       }
404e26ad7d8SSatish Balay     }
4054b3ec233SSatish Balay   }
4064b3ec233SSatish Balay   /* End Send-Recv of row_values */
407*606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
408*606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
4094b3ec233SSatish Balay   s_status2 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status1);
4104b3ec233SSatish Balay   ierr      = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);
411*606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
412*606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
413e26ad7d8SSatish Balay 
414e26ad7d8SSatish Balay   /* Restore the indices */
415e26ad7d8SSatish Balay   for ( i=0; i<ismax; i++ ) {
416e26ad7d8SSatish Balay     ierr = ISRestoreIndices(isrow[i], irow+i);CHKERRQ(ierr);
417e26ad7d8SSatish Balay     ierr = ISRestoreIndices(iscol[i], icol+i);CHKERRQ(ierr);
418e26ad7d8SSatish Balay   }
419e26ad7d8SSatish Balay 
420e26ad7d8SSatish Balay   /* Destroy allocated memory */
421*606d414cSSatish Balay   ierr = PetscFree(irow);CHKERRQ(ierr);
422*606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
423*606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
424e26ad7d8SSatish Balay 
4254b3ec233SSatish Balay 
4264b3ec233SSatish Balay   for ( i=0; i<nrqs; ++i ) {
427*606d414cSSatish Balay     ierr = PetscFree(rbuf2[i]);CHKERRQ(ierr);
4284b3ec233SSatish Balay   }
429*606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
430*606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
431*606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
4324b3ec233SSatish Balay 
433e26ad7d8SSatish Balay   for ( i=0; i<nrqr; ++i ) {
434*606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
435e26ad7d8SSatish Balay   }
436e26ad7d8SSatish Balay 
437*606d414cSSatish Balay   ierr = PetscFree(sbuf2);CHKERRQ(ierr);
438*606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
439e26ad7d8SSatish Balay 
440e26ad7d8SSatish Balay   ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
441e26ad7d8SSatish Balay 
442e26ad7d8SSatish Balay   for ( i=0; i<ismax; i++ ) {
443e26ad7d8SSatish Balay     ierr = MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
444e26ad7d8SSatish Balay     ierr = MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
445e26ad7d8SSatish Balay   }
446e26ad7d8SSatish Balay 
447e26ad7d8SSatish Balay   PetscFunctionReturn(0);
448e26ad7d8SSatish Balay }
449e26ad7d8SSatish Balay 
450