xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision ca44d042d6f86ecc01fa7a52a5213a5161f95f53)
1*ca44d042SBarry Smith /*$Id: mmdense.c,v 1.29 2000/04/12 04:22:59 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;
15d9f96c7cSLois Curfman McInnes   IS           tofrom;
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 */
23029af93fSBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,mdn->N,0,1,&tofrom);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 
28be0abb6dSBarry Smith   ierr = VecCreateMPI(mat->comm,mdn->nvec,mdn->N,&gvec);CHKERRQ(ierr);
29e19c7942SLois Curfman McInnes 
300500aeceSBarry Smith   /* Generate the scatter context */
31d9f96c7cSLois Curfman McInnes   ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx);CHKERRQ(ierr);
320500aeceSBarry Smith   PLogObjectParent(mat,mdn->Mvctx);
330500aeceSBarry Smith   PLogObjectParent(mat,mdn->lvec);
34d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,tofrom);
35d9f96c7cSLois Curfman McInnes   PLogObjectParent(mat,gvec);
36e19c7942SLois Curfman McInnes 
37d9f96c7cSLois Curfman McInnes   ierr = ISDestroy(tofrom);CHKERRQ(ierr);
38d9f96c7cSLois Curfman McInnes   ierr = VecDestroy(gvec);CHKERRQ(ierr);
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
400500aeceSBarry Smith }
410500aeceSBarry Smith 
42*ca44d042SBarry Smith EXTERN int MatGetSubMatrices_MPIDense_Local(Mat,int,IS*,IS*,MatReuse,Mat*);
43e26ad7d8SSatish Balay #undef __FUNC__
44b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrices_MPIDense"
456831982aSBarry Smith int MatGetSubMatrices_MPIDense(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat **submat)
46e26ad7d8SSatish Balay {
47e26ad7d8SSatish Balay   Mat_MPIDense  *c = (Mat_MPIDense*)C->data;
48e26ad7d8SSatish Balay   int           nmax,nstages_local,nstages,i,pos,max_no,ierr;
49e26ad7d8SSatish Balay 
50e26ad7d8SSatish Balay   PetscFunctionBegin;
51e26ad7d8SSatish Balay   /* Allocate memory to hold all the submatrices */
52e26ad7d8SSatish Balay   if (scall != MAT_REUSE_MATRIX) {
53e26ad7d8SSatish Balay     *submat = (Mat *)PetscMalloc((ismax+1)*sizeof(Mat));CHKPTRQ(*submat);
54e26ad7d8SSatish Balay   }
55e26ad7d8SSatish Balay   /* Determine the number of stages through which submatrices are done */
56e26ad7d8SSatish Balay   nmax          = 20*1000000 / (c->N * sizeof(int));
57329f5518SBarry Smith   if (!nmax) nmax = 1;
58e26ad7d8SSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
59e26ad7d8SSatish Balay 
60e26ad7d8SSatish Balay   /* Make sure every processor loops through the nstages */
61e26ad7d8SSatish Balay   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr);
620500aeceSBarry Smith 
630500aeceSBarry Smith 
64e26ad7d8SSatish Balay   for (i=0,pos=0; i<nstages; i++) {
65e26ad7d8SSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
66e26ad7d8SSatish Balay     else if (pos == ismax) max_no = 0;
67e26ad7d8SSatish Balay     else                   max_no = ismax-pos;
68e26ad7d8SSatish Balay     ierr = MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);CHKERRQ(ierr);
69e26ad7d8SSatish Balay     pos += max_no;
70e26ad7d8SSatish Balay   }
71e26ad7d8SSatish Balay   PetscFunctionReturn(0);
72e26ad7d8SSatish Balay }
73e26ad7d8SSatish Balay /* -------------------------------------------------------------------------*/
74e26ad7d8SSatish Balay #undef __FUNC__
75b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrices_MPIDense_Local"
767b2a1423SBarry Smith int MatGetSubMatrices_MPIDense_Local(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat *submats)
77e26ad7d8SSatish Balay {
78e26ad7d8SSatish Balay   Mat_MPIDense  *c = (Mat_MPIDense*)C->data;
79e26ad7d8SSatish Balay   Mat            A = c->A;
80e26ad7d8SSatish Balay   Mat_SeqDense  *a = (Mat_SeqDense*)A->data,*mat;
814b3ec233SSatish Balay   int           N = c->N,rstart = c->rstart,count;
82e26ad7d8SSatish Balay   int           **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size;
834b3ec233SSatish Balay   int           **sbuf1,rank,m,i,j,k,l,ct1,ierr,**rbuf1,row,proc;
844b3ec233SSatish Balay   int           nrqs,msz,**ptr,index,*ctr,*pa,*tmp,bsz,nrqr;
854b3ec233SSatish Balay   int           is_no,jmax,*irow_i,**rmap,*rmap_i;
864b3ec233SSatish Balay   int           len,ctr_j,*sbuf1_j,*rbuf1_i;
874b3ec233SSatish Balay   int           tag0,tag1;
884b3ec233SSatish Balay   MPI_Request   *s_waits1,*r_waits1,*s_waits2,*r_waits2;
894b3ec233SSatish Balay   MPI_Status    *r_status1,*r_status2,*s_status1,*s_status2;
90e26ad7d8SSatish Balay   MPI_Comm      comm;
914b3ec233SSatish Balay   Scalar        **rbuf2,**sbuf2;
92e26ad7d8SSatish Balay 
93e26ad7d8SSatish Balay   PetscFunctionBegin;
94e26ad7d8SSatish Balay   comm   = C->comm;
95e26ad7d8SSatish Balay   tag0   = C->tag;
96e26ad7d8SSatish Balay   size   = c->size;
97e26ad7d8SSatish Balay   rank   = c->rank;
98e26ad7d8SSatish Balay   m      = c->M;
99e26ad7d8SSatish Balay 
100e26ad7d8SSatish Balay   /* Get some new tags to keep the communication clean */
101e26ad7d8SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
102e26ad7d8SSatish Balay 
103e26ad7d8SSatish Balay     /* Check if the col indices are sorted */
104e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
105888f2ed8SSatish Balay     ierr = ISSorted(isrow[i],(PetscTruth*)&j);CHKERRQ(ierr);
106e26ad7d8SSatish Balay     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted");
107888f2ed8SSatish Balay     ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr);
108e26ad7d8SSatish Balay     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted");
109e26ad7d8SSatish Balay   }
110e26ad7d8SSatish Balay 
1114b3ec233SSatish Balay   len    =  2*ismax*(sizeof(int *)+sizeof(int)) + (m+1)*sizeof(int);
112e26ad7d8SSatish Balay   irow   = (int **)PetscMalloc(len);CHKPTRQ(irow);
113e26ad7d8SSatish Balay   icol   = irow + ismax;
114e26ad7d8SSatish Balay   nrow   = (int*)(icol + ismax);
115e26ad7d8SSatish Balay   ncol   = nrow + ismax;
116e26ad7d8SSatish Balay   rtable = ncol + ismax;
117e26ad7d8SSatish Balay 
118e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
119e26ad7d8SSatish Balay     ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
120e26ad7d8SSatish Balay     ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
121e26ad7d8SSatish Balay     ierr = ISGetSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
122e26ad7d8SSatish Balay     ierr = ISGetSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
123e26ad7d8SSatish Balay   }
124e26ad7d8SSatish Balay 
125e26ad7d8SSatish Balay   /* Create hash table for the mapping :row -> proc*/
126e26ad7d8SSatish Balay   for (i=0,j=0; i<size; i++) {
127e26ad7d8SSatish Balay     jmax = c->rowners[i+1];
128e26ad7d8SSatish Balay     for (; j<jmax; j++) {
129e26ad7d8SSatish Balay       rtable[j] = i;
130e26ad7d8SSatish Balay     }
131e26ad7d8SSatish Balay   }
132e26ad7d8SSatish Balay 
133e26ad7d8SSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
134e26ad7d8SSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
135e26ad7d8SSatish Balay   w1     = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1); /* mesg size */
136e26ad7d8SSatish Balay   w2     = w1 + size;      /* if w2[i] marked, then a message to proc i*/
137e26ad7d8SSatish Balay   w3     = w2 + size;      /* no of IS that needs to be sent to proc i */
138e26ad7d8SSatish Balay   w4     = w3 + size;      /* temp work space used in determining w1, w2, w3 */
139549d3d68SSatish Balay   ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/
140e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
141549d3d68SSatish Balay     ierr   = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/
142e26ad7d8SSatish Balay     jmax   = nrow[i];
143e26ad7d8SSatish Balay     irow_i = irow[i];
144e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {
145e26ad7d8SSatish Balay       row  = irow_i[j];
146e26ad7d8SSatish Balay       proc = rtable[row];
147e26ad7d8SSatish Balay       w4[proc]++;
148e26ad7d8SSatish Balay     }
149e26ad7d8SSatish Balay     for (j=0; j<size; j++) {
150e26ad7d8SSatish Balay       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
151e26ad7d8SSatish Balay     }
152e26ad7d8SSatish Balay   }
153e26ad7d8SSatish Balay 
154e26ad7d8SSatish Balay   nrqs     = 0;              /* no of outgoing messages */
155e26ad7d8SSatish Balay   msz      = 0;              /* total mesg length (for all procs) */
156e26ad7d8SSatish Balay   w1[rank] = 0;              /* no mesg sent to self */
157e26ad7d8SSatish Balay   w3[rank] = 0;
158e26ad7d8SSatish Balay   for (i=0; i<size; i++) {
159e26ad7d8SSatish Balay     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
160e26ad7d8SSatish Balay   }
161e26ad7d8SSatish Balay   pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); /*(proc -array)*/
162e26ad7d8SSatish Balay   for (i=0,j=0; i<size; i++) {
163e26ad7d8SSatish Balay     if (w1[i]) { pa[j] = i; j++; }
164e26ad7d8SSatish Balay   }
165e26ad7d8SSatish Balay 
166e26ad7d8SSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
167e26ad7d8SSatish Balay   for (i=0; i<nrqs; i++) {
168e26ad7d8SSatish Balay     j     = pa[i];
169e26ad7d8SSatish Balay     w1[j] += w2[j] + 2* w3[j];
170e26ad7d8SSatish Balay     msz   += w1[j];
171e26ad7d8SSatish Balay   }
172e26ad7d8SSatish Balay   /* Do a global reduction to determine how many messages to expect*/
173e26ad7d8SSatish Balay   {
1746831982aSBarry Smith     int *rw1;
175e26ad7d8SSatish Balay     rw1   = (int *)PetscMalloc(2*size*sizeof(int));CHKPTRQ(rw1);
1766831982aSBarry Smith     ierr  = MPI_Allreduce(w1,rw1,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr);
177e26ad7d8SSatish Balay     bsz   = rw1[rank];
1786831982aSBarry Smith     nrqr  = rw1[size+rank];
179606d414cSSatish Balay     ierr  = PetscFree(rw1);CHKERRQ(ierr);
180e26ad7d8SSatish Balay   }
181e26ad7d8SSatish Balay 
182e26ad7d8SSatish Balay   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
183e26ad7d8SSatish Balay   len      = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
184e26ad7d8SSatish Balay   rbuf1    = (int**)PetscMalloc(len);CHKPTRQ(rbuf1);
185e26ad7d8SSatish Balay   rbuf1[0] = (int*)(rbuf1 + nrqr);
186e26ad7d8SSatish Balay   for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
187e26ad7d8SSatish Balay 
188e26ad7d8SSatish Balay   /* Post the receives */
189e26ad7d8SSatish Balay   r_waits1 = (MPI_Request*)PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(r_waits1);
190e26ad7d8SSatish Balay   for (i=0; i<nrqr; ++i) {
191e26ad7d8SSatish Balay     ierr = MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr);
192e26ad7d8SSatish Balay   }
193e26ad7d8SSatish Balay 
194e26ad7d8SSatish Balay   /* Allocate Memory for outgoing messages */
195e26ad7d8SSatish Balay   len   = 2*size*sizeof(int*) + 2*msz*sizeof(int)+ size*sizeof(int);
196e26ad7d8SSatish Balay   sbuf1 = (int **)PetscMalloc(len);CHKPTRQ(sbuf1);
197e26ad7d8SSatish Balay   ptr   = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
198549d3d68SSatish Balay   ierr  = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr);
199e26ad7d8SSatish Balay   /* allocate memory for outgoing data + buf to receive the first reply */
200e26ad7d8SSatish Balay   tmp   = (int*)(ptr + size);
201e26ad7d8SSatish Balay   ctr   = tmp + 2*msz;
202e26ad7d8SSatish Balay 
203e26ad7d8SSatish Balay   {
204e26ad7d8SSatish Balay     int *iptr = tmp,ict = 0;
205e26ad7d8SSatish Balay     for (i=0; i<nrqs; i++) {
206e26ad7d8SSatish Balay       j         = pa[i];
207e26ad7d8SSatish Balay       iptr     += ict;
208e26ad7d8SSatish Balay       sbuf1[j]  = iptr;
209e26ad7d8SSatish Balay       ict       = w1[j];
210e26ad7d8SSatish Balay     }
211e26ad7d8SSatish Balay   }
212e26ad7d8SSatish Balay 
213e26ad7d8SSatish Balay   /* Form the outgoing messages */
214e26ad7d8SSatish Balay   /* Initialize the header space */
215e26ad7d8SSatish Balay   for (i=0; i<nrqs; i++) {
216e26ad7d8SSatish Balay     j           = pa[i];
217e26ad7d8SSatish Balay     sbuf1[j][0] = 0;
218549d3d68SSatish Balay     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr);
219e26ad7d8SSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
220e26ad7d8SSatish Balay   }
221e26ad7d8SSatish Balay 
222e26ad7d8SSatish Balay   /* Parse the isrow and copy data into outbuf */
223e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
224549d3d68SSatish Balay     ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr);
225e26ad7d8SSatish Balay     irow_i = irow[i];
226e26ad7d8SSatish Balay     jmax   = nrow[i];
227e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
228e26ad7d8SSatish Balay       row  = irow_i[j];
229e26ad7d8SSatish Balay       proc = rtable[row];
230e26ad7d8SSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
231e26ad7d8SSatish Balay         ctr[proc]++;
232e26ad7d8SSatish Balay         *ptr[proc] = row;
233e26ad7d8SSatish Balay         ptr[proc]++;
234e26ad7d8SSatish Balay       }
235e26ad7d8SSatish Balay     }
236e26ad7d8SSatish Balay     /* Update the headers for the current IS */
237e26ad7d8SSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
238e26ad7d8SSatish Balay       if ((ctr_j = ctr[j])) {
239e26ad7d8SSatish Balay         sbuf1_j        = sbuf1[j];
240e26ad7d8SSatish Balay         k              = ++sbuf1_j[0];
241e26ad7d8SSatish Balay         sbuf1_j[2*k]   = ctr_j;
242e26ad7d8SSatish Balay         sbuf1_j[2*k-1] = i;
243e26ad7d8SSatish Balay       }
244e26ad7d8SSatish Balay     }
245e26ad7d8SSatish Balay   }
246e26ad7d8SSatish Balay 
247e26ad7d8SSatish Balay   /*  Now  post the sends */
248e26ad7d8SSatish Balay   s_waits1 = (MPI_Request*)PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(s_waits1);
249e26ad7d8SSatish Balay   for (i=0; i<nrqs; ++i) {
250e26ad7d8SSatish Balay     j = pa[i];
251e26ad7d8SSatish Balay     ierr = MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
252e26ad7d8SSatish Balay   }
253e26ad7d8SSatish Balay 
2544b3ec233SSatish Balay   /* Post recieves to capture the row_data from other procs */
255e26ad7d8SSatish Balay   r_waits2 = (MPI_Request*)PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits2);
2564b3ec233SSatish Balay   rbuf2    = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar *));CHKPTRQ(rbuf2);
2574b3ec233SSatish Balay   for (i=0; i<nrqs; i++) {
258e26ad7d8SSatish Balay     j        = pa[i];
2594b3ec233SSatish Balay     count    = (w1[j] - (2*sbuf1[j][0] + 1))*N;
2604b3ec233SSatish Balay     rbuf2[i] = (Scalar *)PetscMalloc((count+1)*sizeof(Scalar));CHKPTRQ(rbuf2[i]);
2614b3ec233SSatish Balay     ierr     = MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
262e26ad7d8SSatish Balay   }
263e26ad7d8SSatish Balay 
2644b3ec233SSatish Balay   /* Receive messages(row_nos) and then, pack and send off the rowvalues
2654b3ec233SSatish Balay      to the correct processors */
266e26ad7d8SSatish Balay 
267e26ad7d8SSatish Balay   s_waits2  = (MPI_Request*)PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits2);
268e26ad7d8SSatish Balay   r_status1 = (MPI_Status*)PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(r_status1);
2694b3ec233SSatish Balay   sbuf2       = (Scalar**)PetscMalloc((nrqr+1)*sizeof(Scalar*));CHKPTRQ(sbuf2);
270e26ad7d8SSatish Balay 
271e26ad7d8SSatish Balay   {
2724b3ec233SSatish Balay     Scalar *sbuf2_i,*v_start;
2734b3ec233SSatish Balay     int    s_proc;
274e26ad7d8SSatish Balay     for (i=0; i<nrqr; ++i) {
275e26ad7d8SSatish Balay       ierr = MPI_Waitany(nrqr,r_waits1,&index,r_status1+i);CHKERRQ(ierr);
2764b3ec233SSatish Balay       s_proc          = r_status1[i].MPI_SOURCE; /* send processor */
2774b3ec233SSatish Balay       rbuf1_i         = rbuf1[index]; /* Actual message from s_proc */
2784b3ec233SSatish Balay       /* no of rows = end - start; since start is array index[], 0index, whel end
2794b3ec233SSatish Balay          is length of the buffer - which is 1index */
280e26ad7d8SSatish Balay       start           = 2*rbuf1_i[0] + 1;
281e26ad7d8SSatish Balay       MPI_Get_count(r_status1+i,MPI_INT,&end);
2824b3ec233SSatish Balay       /* allocate memory sufficinet to hold all the row values */
2834b3ec233SSatish Balay       sbuf2[index] = (Scalar *)PetscMalloc((end-start)*N*sizeof(Scalar));CHKPTRQ(sbuf2[index]);
284e26ad7d8SSatish Balay       sbuf2_i      = sbuf2[index];
2854b3ec233SSatish Balay       /* Now pack the data */
286e26ad7d8SSatish Balay       for (j=start; j<end; j++) {
2874b3ec233SSatish Balay         row = rbuf1_i[j] - rstart;
2884b3ec233SSatish Balay         v_start = a->v + row;
2894b3ec233SSatish Balay         for (k=0; k<N; k++) {
2904b3ec233SSatish Balay           sbuf2_i[0] = v_start[0];
2914b3ec233SSatish Balay           sbuf2_i++; v_start+=a->m;
292e26ad7d8SSatish Balay         }
293e26ad7d8SSatish Balay       }
2944b3ec233SSatish Balay       /* Now send off the data */
2954b3ec233SSatish Balay       ierr = MPI_Isend(sbuf2[index],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);CHKERRQ(ierr);
2964b3ec233SSatish Balay     }
2974b3ec233SSatish Balay   }
2984b3ec233SSatish Balay   /* End Send-Recv of IS + row_numbers */
299606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
300606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
301e26ad7d8SSatish Balay   s_status1 = (MPI_Status*)PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(s_status1);
302e26ad7d8SSatish Balay   ierr      = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);
303606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
304606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
305e26ad7d8SSatish Balay 
3064b3ec233SSatish Balay   /* Create the submatrices */
3074b3ec233SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
308e26ad7d8SSatish Balay     for (i=0; i<ismax; i++) {
3094b3ec233SSatish Balay       mat = (Mat_SeqDense *)(submats[i]->data);
3104b3ec233SSatish Balay       if ((mat->m != nrow[i]) || (mat->n != ncol[i])) {
3114b3ec233SSatish Balay         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
312e26ad7d8SSatish Balay       }
313549d3d68SSatish Balay       ierr = PetscMemzero(mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr);
3144b3ec233SSatish Balay       submats[i]->factor = C->factor;
315e26ad7d8SSatish Balay     }
3164b3ec233SSatish Balay   } else {
317e26ad7d8SSatish Balay     for (i=0; i<ismax; i++) {
3184b3ec233SSatish Balay       ierr = MatCreateSeqDense(PETSC_COMM_SELF,nrow[i],ncol[i],PETSC_NULL,submats+i);CHKERRQ(ierr);
3194b3ec233SSatish Balay     }
3204b3ec233SSatish Balay   }
3214b3ec233SSatish Balay 
3224b3ec233SSatish Balay   /* Assemble the matrices */
3234b3ec233SSatish Balay   {
3244b3ec233SSatish Balay     int    col;
3254b3ec233SSatish Balay     Scalar *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];
3424b3ec233SSatish Balay             imat_vi[k*m] = mat_vi[col*a->m];
343e26ad7d8SSatish Balay           }
3444b3ec233SSatish Balay         }
345e26ad7d8SSatish Balay       }
346e26ad7d8SSatish Balay     }
347e26ad7d8SSatish Balay   }
348e26ad7d8SSatish Balay 
3494b3ec233SSatish Balay   /* Create row map. This maps c->row to submat->row for each submat*/
3504b3ec233SSatish Balay   /* this is a very expensive operation wrt memory usage */
351e26ad7d8SSatish Balay   len     = (1+ismax)*sizeof(int*)+ ismax*c->M*sizeof(int);
352e26ad7d8SSatish Balay   rmap    = (int **)PetscMalloc(len);CHKPTRQ(rmap);
353e26ad7d8SSatish Balay   rmap[0] = (int *)(rmap + ismax);
354549d3d68SSatish Balay   ierr    = PetscMemzero(rmap[0],ismax*c->M*sizeof(int));CHKERRQ(ierr);
355e26ad7d8SSatish Balay   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + c->M;}
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 */
366e26ad7d8SSatish Balay 
3674b3ec233SSatish Balay   r_status2 = (MPI_Status*)PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status2);
3684b3ec233SSatish Balay 
3694b3ec233SSatish Balay   {
3704b3ec233SSatish Balay     int    is_max,tmp1,col,*sbuf1_i,is_sz;
3714b3ec233SSatish Balay     Scalar *rbuf2_i,*imat_v,*imat_vi;
3724b3ec233SSatish Balay 
3734b3ec233SSatish Balay     for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */
3744b3ec233SSatish Balay       ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr);
3754b3ec233SSatish Balay       /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
3764b3ec233SSatish Balay       sbuf1_i = sbuf1[pa[i]];
3774b3ec233SSatish Balay       is_max  = sbuf1_i[0];
3784b3ec233SSatish Balay       ct1     = 2*is_max+1;
379e26ad7d8SSatish Balay       rbuf2_i = rbuf2[i];
3804b3ec233SSatish Balay       for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */
381e26ad7d8SSatish Balay         is_no     = sbuf1_i[2*j-1];
3824b3ec233SSatish Balay         is_sz     = sbuf1_i[2*j];
383e26ad7d8SSatish Balay         mat       = (Mat_SeqDense*)submats[is_no]->data;
3844b3ec233SSatish Balay         imat_v    = mat->v;
3854b3ec233SSatish Balay         rmap_i    = rmap[is_no];
3864b3ec233SSatish Balay         m         = nrow[is_no];
3874b3ec233SSatish Balay         for (k=0; k<is_sz; k++,rbuf2_i+=N) {  /* For each row */
3884b3ec233SSatish Balay           row      = sbuf1_i[ct1]; ct1++;
389e26ad7d8SSatish Balay           row      = rmap_i[row];
3904b3ec233SSatish Balay           imat_vi  = imat_v + row;
3914b3ec233SSatish Balay           for (l=0; l<ncol[is_no]; l++) { /* For each col */
3924b3ec233SSatish Balay             col = icol[is_no][l];
3934b3ec233SSatish Balay             imat_vi[l*m] = rbuf2_i[col];
394e26ad7d8SSatish Balay           }
395e26ad7d8SSatish Balay         }
396e26ad7d8SSatish Balay       }
397e26ad7d8SSatish Balay     }
3984b3ec233SSatish Balay   }
3994b3ec233SSatish Balay   /* End Send-Recv of row_values */
400606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
401606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
4024b3ec233SSatish Balay   s_status2 = (MPI_Status*)PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status1);
4034b3ec233SSatish Balay   ierr      = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);
404606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
405606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
406e26ad7d8SSatish Balay 
407e26ad7d8SSatish Balay   /* Restore the indices */
408e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
409e26ad7d8SSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
410e26ad7d8SSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
411e26ad7d8SSatish Balay   }
412e26ad7d8SSatish Balay 
413e26ad7d8SSatish Balay   /* Destroy allocated memory */
414606d414cSSatish Balay   ierr = PetscFree(irow);CHKERRQ(ierr);
415606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
416606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
417e26ad7d8SSatish Balay 
4184b3ec233SSatish Balay 
4194b3ec233SSatish Balay   for (i=0; i<nrqs; ++i) {
420606d414cSSatish Balay     ierr = PetscFree(rbuf2[i]);CHKERRQ(ierr);
4214b3ec233SSatish Balay   }
422606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
423606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
424606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
4254b3ec233SSatish Balay 
426e26ad7d8SSatish Balay   for (i=0; i<nrqr; ++i) {
427606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
428e26ad7d8SSatish Balay   }
429e26ad7d8SSatish Balay 
430606d414cSSatish Balay   ierr = PetscFree(sbuf2);CHKERRQ(ierr);
431606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
432e26ad7d8SSatish Balay 
433e26ad7d8SSatish Balay   ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
434e26ad7d8SSatish Balay 
435e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
436e26ad7d8SSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
437e26ad7d8SSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
438e26ad7d8SSatish Balay   }
439e26ad7d8SSatish Balay 
440e26ad7d8SSatish Balay   PetscFunctionReturn(0);
441e26ad7d8SSatish Balay }
442e26ad7d8SSatish Balay 
443