xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision b9b9770376f6bdad54cc53bb7a8b79bbb53b7e80)
1*b9b97703SBarry Smith /*$Id: mmdense.c,v 1.30 2000/05/10 16:40:34 bsmith Exp bsmith $*/
20500aeceSBarry Smith 
30500aeceSBarry Smith /*
40500aeceSBarry Smith    Support for the parallel dense matrix vector multiply
50500aeceSBarry Smith */
670f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
7f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
80500aeceSBarry Smith 
95615d1e5SSatish Balay #undef __FUNC__
10b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetUpMultiply_MPIDense"
110500aeceSBarry Smith int MatSetUpMultiply_MPIDense(Mat mat)
120500aeceSBarry Smith {
130500aeceSBarry Smith   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
146f0b7f37SSatish Balay   int          ierr;
15*b9b97703SBarry Smith   IS           from,to;
160500aeceSBarry Smith   Vec          gvec;
170500aeceSBarry Smith 
183a40ed3dSBarry Smith   PetscFunctionBegin;
190500aeceSBarry Smith   /* Create local vector that is used to scatter into */
20029af93fSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,mdn->N,&mdn->lvec);CHKERRQ(ierr);
210500aeceSBarry Smith 
220500aeceSBarry Smith   /* Create temporary index set for building scatter gather */
23*b9b97703SBarry Smith   ierr = ISCreateStride(mat->comm,mdn->N,0,1,&from);CHKERRQ(ierr);
24*b9b97703SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,mdn->N,0,1,&to);CHKERRQ(ierr);
250500aeceSBarry Smith 
260500aeceSBarry Smith   /* Create temporary global vector to generate scatter context */
276f0b7f37SSatish Balay   /* n    = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank]; */
28426b7fd2SBarry Smith 
29be0abb6dSBarry Smith   ierr = VecCreateMPI(mat->comm,mdn->nvec,mdn->N,&gvec);CHKERRQ(ierr);
30e19c7942SLois Curfman McInnes 
310500aeceSBarry Smith   /* Generate the scatter context */
32*b9b97703SBarry Smith   ierr = VecScatterCreate(gvec,from,mdn->lvec,to,&mdn->Mvctx);CHKERRQ(ierr);
330500aeceSBarry Smith   PLogObjectParent(mat,mdn->Mvctx);
340500aeceSBarry Smith   PLogObjectParent(mat,mdn->lvec);
35*b9b97703SBarry Smith   PLogObjectParent(mat,from);
36*b9b97703SBarry Smith   PLogObjectParent(mat,to);
37d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,gvec);
38e19c7942SLois Curfman McInnes 
39*b9b97703SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
40*b9b97703SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
41d9f96c7cSLois Curfman McInnes   ierr = VecDestroy(gvec);CHKERRQ(ierr);
423a40ed3dSBarry Smith   PetscFunctionReturn(0);
430500aeceSBarry Smith }
440500aeceSBarry Smith 
45ca44d042SBarry Smith EXTERN int MatGetSubMatrices_MPIDense_Local(Mat,int,IS*,IS*,MatReuse,Mat*);
46e26ad7d8SSatish Balay #undef __FUNC__
47b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrices_MPIDense"
486831982aSBarry Smith int MatGetSubMatrices_MPIDense(Mat C,int ismax,IS *isrow,IS *iscol,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));
60329f5518SBarry Smith   if (!nmax) 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__
78b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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);
124*b9b97703SBarry Smith     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
125*b9b97703SBarry Smith     ierr = ISGetLocalSize(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   {
1776831982aSBarry Smith     int *rw1;
178e26ad7d8SSatish Balay     rw1   = (int *)PetscMalloc(2*size*sizeof(int));CHKPTRQ(rw1);
1796831982aSBarry Smith     ierr  = MPI_Allreduce(w1,rw1,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr);
180e26ad7d8SSatish Balay     bsz   = rw1[rank];
1816831982aSBarry Smith     nrqr  = rw1[size+rank];
182606d414cSSatish Balay     ierr  = PetscFree(rw1);CHKERRQ(ierr);
183e26ad7d8SSatish Balay   }
184e26ad7d8SSatish Balay 
185e26ad7d8SSatish Balay   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
186e26ad7d8SSatish Balay   len      = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
187e26ad7d8SSatish Balay   rbuf1    = (int**)PetscMalloc(len);CHKPTRQ(rbuf1);
188e26ad7d8SSatish Balay   rbuf1[0] = (int*)(rbuf1 + nrqr);
189e26ad7d8SSatish Balay   for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
190e26ad7d8SSatish Balay 
191e26ad7d8SSatish Balay   /* Post the receives */
192e26ad7d8SSatish Balay   r_waits1 = (MPI_Request*)PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(r_waits1);
193e26ad7d8SSatish Balay   for (i=0; i<nrqr; ++i) {
194e26ad7d8SSatish Balay     ierr = MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr);
195e26ad7d8SSatish Balay   }
196e26ad7d8SSatish Balay 
197e26ad7d8SSatish Balay   /* Allocate Memory for outgoing messages */
198e26ad7d8SSatish Balay   len   = 2*size*sizeof(int*) + 2*msz*sizeof(int)+ size*sizeof(int);
199e26ad7d8SSatish Balay   sbuf1 = (int **)PetscMalloc(len);CHKPTRQ(sbuf1);
200e26ad7d8SSatish Balay   ptr   = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
201549d3d68SSatish Balay   ierr  = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr);
202e26ad7d8SSatish Balay   /* allocate memory for outgoing data + buf to receive the first reply */
203e26ad7d8SSatish Balay   tmp   = (int*)(ptr + size);
204e26ad7d8SSatish Balay   ctr   = tmp + 2*msz;
205e26ad7d8SSatish Balay 
206e26ad7d8SSatish Balay   {
207e26ad7d8SSatish Balay     int *iptr = tmp,ict = 0;
208e26ad7d8SSatish Balay     for (i=0; i<nrqs; i++) {
209e26ad7d8SSatish Balay       j         = pa[i];
210e26ad7d8SSatish Balay       iptr     += ict;
211e26ad7d8SSatish Balay       sbuf1[j]  = iptr;
212e26ad7d8SSatish Balay       ict       = w1[j];
213e26ad7d8SSatish Balay     }
214e26ad7d8SSatish Balay   }
215e26ad7d8SSatish Balay 
216e26ad7d8SSatish Balay   /* Form the outgoing messages */
217e26ad7d8SSatish Balay   /* Initialize the header space */
218e26ad7d8SSatish Balay   for (i=0; i<nrqs; i++) {
219e26ad7d8SSatish Balay     j           = pa[i];
220e26ad7d8SSatish Balay     sbuf1[j][0] = 0;
221549d3d68SSatish Balay     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr);
222e26ad7d8SSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
223e26ad7d8SSatish Balay   }
224e26ad7d8SSatish Balay 
225e26ad7d8SSatish Balay   /* Parse the isrow and copy data into outbuf */
226e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
227549d3d68SSatish Balay     ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr);
228e26ad7d8SSatish Balay     irow_i = irow[i];
229e26ad7d8SSatish Balay     jmax   = nrow[i];
230e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
231e26ad7d8SSatish Balay       row  = irow_i[j];
232e26ad7d8SSatish Balay       proc = rtable[row];
233e26ad7d8SSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
234e26ad7d8SSatish Balay         ctr[proc]++;
235e26ad7d8SSatish Balay         *ptr[proc] = row;
236e26ad7d8SSatish Balay         ptr[proc]++;
237e26ad7d8SSatish Balay       }
238e26ad7d8SSatish Balay     }
239e26ad7d8SSatish Balay     /* Update the headers for the current IS */
240e26ad7d8SSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
241e26ad7d8SSatish Balay       if ((ctr_j = ctr[j])) {
242e26ad7d8SSatish Balay         sbuf1_j        = sbuf1[j];
243e26ad7d8SSatish Balay         k              = ++sbuf1_j[0];
244e26ad7d8SSatish Balay         sbuf1_j[2*k]   = ctr_j;
245e26ad7d8SSatish Balay         sbuf1_j[2*k-1] = i;
246e26ad7d8SSatish Balay       }
247e26ad7d8SSatish Balay     }
248e26ad7d8SSatish Balay   }
249e26ad7d8SSatish Balay 
250e26ad7d8SSatish Balay   /*  Now  post the sends */
251e26ad7d8SSatish Balay   s_waits1 = (MPI_Request*)PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(s_waits1);
252e26ad7d8SSatish Balay   for (i=0; i<nrqs; ++i) {
253e26ad7d8SSatish Balay     j = pa[i];
254e26ad7d8SSatish Balay     ierr = MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
255e26ad7d8SSatish Balay   }
256e26ad7d8SSatish Balay 
2574b3ec233SSatish Balay   /* Post recieves to capture the row_data from other procs */
258e26ad7d8SSatish Balay   r_waits2 = (MPI_Request*)PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits2);
2594b3ec233SSatish Balay   rbuf2    = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar *));CHKPTRQ(rbuf2);
2604b3ec233SSatish Balay   for (i=0; i<nrqs; i++) {
261e26ad7d8SSatish Balay     j        = pa[i];
2624b3ec233SSatish Balay     count    = (w1[j] - (2*sbuf1[j][0] + 1))*N;
2634b3ec233SSatish Balay     rbuf2[i] = (Scalar *)PetscMalloc((count+1)*sizeof(Scalar));CHKPTRQ(rbuf2[i]);
2644b3ec233SSatish Balay     ierr     = MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
265e26ad7d8SSatish Balay   }
266e26ad7d8SSatish Balay 
2674b3ec233SSatish Balay   /* Receive messages(row_nos) and then, pack and send off the rowvalues
2684b3ec233SSatish Balay      to the correct processors */
269e26ad7d8SSatish Balay 
270e26ad7d8SSatish Balay   s_waits2  = (MPI_Request*)PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits2);
271e26ad7d8SSatish Balay   r_status1 = (MPI_Status*)PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(r_status1);
2724b3ec233SSatish Balay   sbuf2       = (Scalar**)PetscMalloc((nrqr+1)*sizeof(Scalar*));CHKPTRQ(sbuf2);
273e26ad7d8SSatish Balay 
274e26ad7d8SSatish Balay   {
2754b3ec233SSatish Balay     Scalar *sbuf2_i,*v_start;
2764b3ec233SSatish Balay     int    s_proc;
277e26ad7d8SSatish Balay     for (i=0; i<nrqr; ++i) {
278e26ad7d8SSatish Balay       ierr = MPI_Waitany(nrqr,r_waits1,&index,r_status1+i);CHKERRQ(ierr);
2794b3ec233SSatish Balay       s_proc          = r_status1[i].MPI_SOURCE; /* send processor */
2804b3ec233SSatish Balay       rbuf1_i         = rbuf1[index]; /* Actual message from s_proc */
2814b3ec233SSatish Balay       /* no of rows = end - start; since start is array index[], 0index, whel end
2824b3ec233SSatish Balay          is length of the buffer - which is 1index */
283e26ad7d8SSatish Balay       start           = 2*rbuf1_i[0] + 1;
284e26ad7d8SSatish Balay       MPI_Get_count(r_status1+i,MPI_INT,&end);
2854b3ec233SSatish Balay       /* allocate memory sufficinet to hold all the row values */
2864b3ec233SSatish Balay       sbuf2[index] = (Scalar *)PetscMalloc((end-start)*N*sizeof(Scalar));CHKPTRQ(sbuf2[index]);
287e26ad7d8SSatish Balay       sbuf2_i      = sbuf2[index];
2884b3ec233SSatish Balay       /* Now pack the data */
289e26ad7d8SSatish Balay       for (j=start; j<end; j++) {
2904b3ec233SSatish Balay         row = rbuf1_i[j] - rstart;
2914b3ec233SSatish Balay         v_start = a->v + row;
2924b3ec233SSatish Balay         for (k=0; k<N; k++) {
2934b3ec233SSatish Balay           sbuf2_i[0] = v_start[0];
2944b3ec233SSatish Balay           sbuf2_i++; v_start+=a->m;
295e26ad7d8SSatish Balay         }
296e26ad7d8SSatish Balay       }
2974b3ec233SSatish Balay       /* Now send off the data */
2984b3ec233SSatish Balay       ierr = MPI_Isend(sbuf2[index],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);CHKERRQ(ierr);
2994b3ec233SSatish Balay     }
3004b3ec233SSatish Balay   }
3014b3ec233SSatish Balay   /* End Send-Recv of IS + row_numbers */
302606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
303606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
304e26ad7d8SSatish Balay   s_status1 = (MPI_Status*)PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(s_status1);
305e26ad7d8SSatish Balay   ierr      = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);
306606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
307606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
308e26ad7d8SSatish Balay 
3094b3ec233SSatish Balay   /* Create the submatrices */
3104b3ec233SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
311e26ad7d8SSatish Balay     for (i=0; i<ismax; i++) {
3124b3ec233SSatish Balay       mat = (Mat_SeqDense *)(submats[i]->data);
3134b3ec233SSatish Balay       if ((mat->m != nrow[i]) || (mat->n != ncol[i])) {
3144b3ec233SSatish Balay         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
315e26ad7d8SSatish Balay       }
316549d3d68SSatish Balay       ierr = PetscMemzero(mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr);
3174b3ec233SSatish Balay       submats[i]->factor = C->factor;
318e26ad7d8SSatish Balay     }
3194b3ec233SSatish Balay   } else {
320e26ad7d8SSatish Balay     for (i=0; i<ismax; i++) {
3214b3ec233SSatish Balay       ierr = MatCreateSeqDense(PETSC_COMM_SELF,nrow[i],ncol[i],PETSC_NULL,submats+i);CHKERRQ(ierr);
3224b3ec233SSatish Balay     }
3234b3ec233SSatish Balay   }
3244b3ec233SSatish Balay 
3254b3ec233SSatish Balay   /* Assemble the matrices */
3264b3ec233SSatish Balay   {
3274b3ec233SSatish Balay     int    col;
3284b3ec233SSatish Balay     Scalar *imat_v,*mat_v,*imat_vi,*mat_vi;
3294b3ec233SSatish Balay 
3304b3ec233SSatish Balay     for (i=0; i<ismax; i++) {
3314b3ec233SSatish Balay       mat       = (Mat_SeqDense*)submats[i]->data;
3324b3ec233SSatish Balay       mat_v     = a->v;
3334b3ec233SSatish Balay       imat_v    = mat->v;
334e26ad7d8SSatish Balay       irow_i    = irow[i];
3354b3ec233SSatish Balay       m         = nrow[i];
3364b3ec233SSatish Balay       for (j=0; j<m; j++) {
337e26ad7d8SSatish Balay         row      = irow_i[j] ;
338e26ad7d8SSatish Balay         proc     = rtable[row];
339e26ad7d8SSatish Balay         if (proc == rank) {
3404b3ec233SSatish Balay           row      = row - rstart;
3414b3ec233SSatish Balay           mat_vi   = mat_v + row;
3424b3ec233SSatish Balay           imat_vi  = imat_v + j;
3434b3ec233SSatish Balay           for (k=0; k<ncol[i]; k++) {
3444b3ec233SSatish Balay             col = icol[i][k];
3454b3ec233SSatish Balay             imat_vi[k*m] = mat_vi[col*a->m];
346e26ad7d8SSatish Balay           }
3474b3ec233SSatish Balay         }
348e26ad7d8SSatish Balay       }
349e26ad7d8SSatish Balay     }
350e26ad7d8SSatish Balay   }
351e26ad7d8SSatish Balay 
3524b3ec233SSatish Balay   /* Create row map. This maps c->row to submat->row for each submat*/
3534b3ec233SSatish Balay   /* this is a very expensive operation wrt memory usage */
354e26ad7d8SSatish Balay   len     = (1+ismax)*sizeof(int*)+ ismax*c->M*sizeof(int);
355e26ad7d8SSatish Balay   rmap    = (int **)PetscMalloc(len);CHKPTRQ(rmap);
356e26ad7d8SSatish Balay   rmap[0] = (int *)(rmap + ismax);
357549d3d68SSatish Balay   ierr    = PetscMemzero(rmap[0],ismax*c->M*sizeof(int));CHKERRQ(ierr);
358e26ad7d8SSatish Balay   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + c->M;}
359e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
360e26ad7d8SSatish Balay     rmap_i = rmap[i];
361e26ad7d8SSatish Balay     irow_i = irow[i];
362e26ad7d8SSatish Balay     jmax   = nrow[i];
363e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {
364e26ad7d8SSatish Balay       rmap_i[irow_i[j]] = j;
365e26ad7d8SSatish Balay     }
366e26ad7d8SSatish Balay   }
367e26ad7d8SSatish Balay 
3684b3ec233SSatish Balay   /* Now Receive the row_values and assemble the rest of the matrix */
369e26ad7d8SSatish Balay 
3704b3ec233SSatish Balay   r_status2 = (MPI_Status*)PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status2);
3714b3ec233SSatish Balay 
3724b3ec233SSatish Balay   {
3734b3ec233SSatish Balay     int    is_max,tmp1,col,*sbuf1_i,is_sz;
3744b3ec233SSatish Balay     Scalar *rbuf2_i,*imat_v,*imat_vi;
3754b3ec233SSatish Balay 
3764b3ec233SSatish Balay     for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */
3774b3ec233SSatish Balay       ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr);
3784b3ec233SSatish Balay       /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
3794b3ec233SSatish Balay       sbuf1_i = sbuf1[pa[i]];
3804b3ec233SSatish Balay       is_max  = sbuf1_i[0];
3814b3ec233SSatish Balay       ct1     = 2*is_max+1;
382e26ad7d8SSatish Balay       rbuf2_i = rbuf2[i];
3834b3ec233SSatish Balay       for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */
384e26ad7d8SSatish Balay         is_no     = sbuf1_i[2*j-1];
3854b3ec233SSatish Balay         is_sz     = sbuf1_i[2*j];
386e26ad7d8SSatish Balay         mat       = (Mat_SeqDense*)submats[is_no]->data;
3874b3ec233SSatish Balay         imat_v    = mat->v;
3884b3ec233SSatish Balay         rmap_i    = rmap[is_no];
3894b3ec233SSatish Balay         m         = nrow[is_no];
3904b3ec233SSatish Balay         for (k=0; k<is_sz; k++,rbuf2_i+=N) {  /* For each row */
3914b3ec233SSatish Balay           row      = sbuf1_i[ct1]; ct1++;
392e26ad7d8SSatish Balay           row      = rmap_i[row];
3934b3ec233SSatish Balay           imat_vi  = imat_v + row;
3944b3ec233SSatish Balay           for (l=0; l<ncol[is_no]; l++) { /* For each col */
3954b3ec233SSatish Balay             col = icol[is_no][l];
3964b3ec233SSatish Balay             imat_vi[l*m] = rbuf2_i[col];
397e26ad7d8SSatish Balay           }
398e26ad7d8SSatish Balay         }
399e26ad7d8SSatish Balay       }
400e26ad7d8SSatish Balay     }
4014b3ec233SSatish Balay   }
4024b3ec233SSatish Balay   /* End Send-Recv of row_values */
403606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
404606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
4054b3ec233SSatish Balay   s_status2 = (MPI_Status*)PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status1);
4064b3ec233SSatish Balay   ierr      = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);
407606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
408606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
409e26ad7d8SSatish Balay 
410e26ad7d8SSatish Balay   /* Restore the indices */
411e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
412e26ad7d8SSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
413e26ad7d8SSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
414e26ad7d8SSatish Balay   }
415e26ad7d8SSatish Balay 
416e26ad7d8SSatish Balay   /* Destroy allocated memory */
417606d414cSSatish Balay   ierr = PetscFree(irow);CHKERRQ(ierr);
418606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
419606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
420e26ad7d8SSatish Balay 
4214b3ec233SSatish Balay 
4224b3ec233SSatish Balay   for (i=0; i<nrqs; ++i) {
423606d414cSSatish Balay     ierr = PetscFree(rbuf2[i]);CHKERRQ(ierr);
4244b3ec233SSatish Balay   }
425606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
426606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
427606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
4284b3ec233SSatish Balay 
429e26ad7d8SSatish Balay   for (i=0; i<nrqr; ++i) {
430606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
431e26ad7d8SSatish Balay   }
432e26ad7d8SSatish Balay 
433606d414cSSatish Balay   ierr = PetscFree(sbuf2);CHKERRQ(ierr);
434606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
435e26ad7d8SSatish Balay 
436e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
437e26ad7d8SSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
438e26ad7d8SSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
439e26ad7d8SSatish Balay   }
440e26ad7d8SSatish Balay 
441e26ad7d8SSatish Balay   PetscFunctionReturn(0);
442e26ad7d8SSatish Balay }
443e26ad7d8SSatish Balay 
444