xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
10500aeceSBarry Smith /*
20500aeceSBarry Smith    Support for the parallel dense matrix vector multiply
30500aeceSBarry Smith */
470f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
50500aeceSBarry Smith 
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIDense"
8dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIDense(Mat mat)
90500aeceSBarry Smith {
100500aeceSBarry Smith   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
11dfbe8321SBarry Smith   PetscErrorCode ierr;
12b9b97703SBarry Smith   IS           from,to;
130500aeceSBarry Smith   Vec          gvec;
140500aeceSBarry Smith 
153a40ed3dSBarry Smith   PetscFunctionBegin;
160500aeceSBarry Smith   /* Create local vector that is used to scatter into */
17273d9f13SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,mat->N,&mdn->lvec);CHKERRQ(ierr);
180500aeceSBarry Smith 
190500aeceSBarry Smith   /* Create temporary index set for building scatter gather */
20273d9f13SBarry Smith   ierr = ISCreateStride(mat->comm,mat->N,0,1,&from);CHKERRQ(ierr);
21273d9f13SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,mat->N,0,1,&to);CHKERRQ(ierr);
220500aeceSBarry Smith 
230500aeceSBarry Smith   /* Create temporary global vector to generate scatter context */
246f0b7f37SSatish Balay   /* n    = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank]; */
25426b7fd2SBarry Smith 
26273d9f13SBarry Smith   ierr = VecCreateMPI(mat->comm,mdn->nvec,mat->N,&gvec);CHKERRQ(ierr);
27e19c7942SLois Curfman McInnes 
280500aeceSBarry Smith   /* Generate the scatter context */
29b9b97703SBarry Smith   ierr = VecScatterCreate(gvec,from,mdn->lvec,to,&mdn->Mvctx);CHKERRQ(ierr);
30b0a32e0cSBarry Smith   PetscLogObjectParent(mat,mdn->Mvctx);
31b0a32e0cSBarry Smith   PetscLogObjectParent(mat,mdn->lvec);
32b0a32e0cSBarry Smith   PetscLogObjectParent(mat,from);
33b0a32e0cSBarry Smith   PetscLogObjectParent(mat,to);
34b0a32e0cSBarry Smith   PetscLogObjectParent(mat,gvec);
35e19c7942SLois Curfman McInnes 
36b9b97703SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
37b9b97703SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
38d9f96c7cSLois Curfman McInnes   ierr = VecDestroy(gvec);CHKERRQ(ierr);
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
400500aeceSBarry Smith }
410500aeceSBarry Smith 
42dfbe8321SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_MPIDense_Local(Mat,int,const IS[],const IS[],MatReuse,Mat*);
434a2ae208SSatish Balay #undef __FUNCT__
444a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIDense"
45dfbe8321SBarry Smith PetscErrorCode MatGetSubMatrices_MPIDense(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
46e26ad7d8SSatish Balay {
47*6849ba73SBarry Smith   PetscErrorCode ierr;
48*6849ba73SBarry Smith   int           nmax,nstages_local,nstages,i,pos,max_no;
49e26ad7d8SSatish Balay 
50e26ad7d8SSatish Balay   PetscFunctionBegin;
51e26ad7d8SSatish Balay   /* Allocate memory to hold all the submatrices */
52e26ad7d8SSatish Balay   if (scall != MAT_REUSE_MATRIX) {
53b0a32e0cSBarry Smith     ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr);
54e26ad7d8SSatish Balay   }
55e26ad7d8SSatish Balay   /* Determine the number of stages through which submatrices are done */
56273d9f13SBarry Smith   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 /* -------------------------------------------------------------------------*/
744a2ae208SSatish Balay #undef __FUNCT__
754a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIDense_Local"
76dfbe8321SBarry Smith PetscErrorCode MatGetSubMatrices_MPIDense_Local(Mat C,int ismax,const IS isrow[],const 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;
81*6849ba73SBarry Smith   PetscErrorCode ierr;
82273d9f13SBarry Smith   int           N = C->N,rstart = c->rstart,count;
83c1dc657dSBarry Smith   int           **irow,**icol,*nrow,*ncol,*w1,*w3,*w4,*rtable,start,end,size;
84*6849ba73SBarry Smith   int           **sbuf1,rank,m,i,j,k,l,ct1,**rbuf1,row,proc;
85999d9058SBarry Smith   int           nrqs,msz,**ptr,idex,*ctr,*pa,*tmp,bsz,nrqr;
864b3ec233SSatish Balay   int           is_no,jmax,*irow_i,**rmap,*rmap_i;
874b3ec233SSatish Balay   int           len,ctr_j,*sbuf1_j,*rbuf1_i;
884b3ec233SSatish Balay   int           tag0,tag1;
894b3ec233SSatish Balay   MPI_Request   *s_waits1,*r_waits1,*s_waits2,*r_waits2;
904b3ec233SSatish Balay   MPI_Status    *r_status1,*r_status2,*s_status1,*s_status2;
91e26ad7d8SSatish Balay   MPI_Comm      comm;
9287828ca2SBarry Smith   PetscScalar   **rbuf2,**sbuf2;
935c5985e7SKris Buschelman   PetscTruth    sorted;
94e26ad7d8SSatish Balay 
95e26ad7d8SSatish Balay   PetscFunctionBegin;
96e26ad7d8SSatish Balay   comm   = C->comm;
97e26ad7d8SSatish Balay   tag0   = C->tag;
98e26ad7d8SSatish Balay   size   = c->size;
99e26ad7d8SSatish Balay   rank   = c->rank;
100273d9f13SBarry Smith   m      = C->M;
101e26ad7d8SSatish Balay 
102e26ad7d8SSatish Balay   /* Get some new tags to keep the communication clean */
103e26ad7d8SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
104e26ad7d8SSatish Balay 
105e26ad7d8SSatish Balay     /* Check if the col indices are sorted */
106e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
1075c5985e7SKris Buschelman     ierr = ISSorted(isrow[i],&sorted);CHKERRQ(ierr);
1085c5985e7SKris Buschelman     if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1095c5985e7SKris Buschelman     ierr = ISSorted(iscol[i],&sorted);CHKERRQ(ierr);
1105c5985e7SKris Buschelman     if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
111e26ad7d8SSatish Balay   }
112e26ad7d8SSatish Balay 
1134b3ec233SSatish Balay   len    =  2*ismax*(sizeof(int*)+sizeof(int)) + (m+1)*sizeof(int);
114b0a32e0cSBarry Smith   ierr   = PetscMalloc(len,&irow);CHKERRQ(ierr);
115e26ad7d8SSatish Balay   icol   = irow + ismax;
116e26ad7d8SSatish Balay   nrow   = (int*)(icol + ismax);
117e26ad7d8SSatish Balay   ncol   = nrow + ismax;
118e26ad7d8SSatish Balay   rtable = ncol + ismax;
119e26ad7d8SSatish Balay 
120e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
121e26ad7d8SSatish Balay     ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
122e26ad7d8SSatish Balay     ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
123b9b97703SBarry Smith     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
124b9b97703SBarry Smith     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
125e26ad7d8SSatish Balay   }
126e26ad7d8SSatish Balay 
127e26ad7d8SSatish Balay   /* Create hash table for the mapping :row -> proc*/
128e26ad7d8SSatish Balay   for (i=0,j=0; i<size; i++) {
129e26ad7d8SSatish Balay     jmax = c->rowners[i+1];
130e26ad7d8SSatish Balay     for (; j<jmax; j++) {
131e26ad7d8SSatish Balay       rtable[j] = i;
132e26ad7d8SSatish Balay     }
133e26ad7d8SSatish Balay   }
134e26ad7d8SSatish Balay 
135e26ad7d8SSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
136e26ad7d8SSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
137b0a32e0cSBarry Smith   ierr   = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr); /* mesg size */
138c1dc657dSBarry Smith   w3     = w1 + 2*size;      /* no of IS that needs to be sent to proc i */
139c1dc657dSBarry Smith   w4     = w3 + size;      /* temp work space used in determining w1,  w3 */
140549d3d68SSatish Balay   ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/
141e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
142549d3d68SSatish Balay     ierr   = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/
143e26ad7d8SSatish Balay     jmax   = nrow[i];
144e26ad7d8SSatish Balay     irow_i = irow[i];
145e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {
146e26ad7d8SSatish Balay       row  = irow_i[j];
147e26ad7d8SSatish Balay       proc = rtable[row];
148e26ad7d8SSatish Balay       w4[proc]++;
149e26ad7d8SSatish Balay     }
150e26ad7d8SSatish Balay     for (j=0; j<size; j++) {
151c1dc657dSBarry Smith       if (w4[j]) { w1[2*j] += w4[j];  w3[j]++;}
152e26ad7d8SSatish Balay     }
153e26ad7d8SSatish Balay   }
154e26ad7d8SSatish Balay 
155e26ad7d8SSatish Balay   nrqs       = 0;              /* no of outgoing messages */
156e26ad7d8SSatish Balay   msz        = 0;              /* total mesg length (for all procs) */
157c1dc657dSBarry Smith   w1[2*rank] = 0;              /* no mesg sent to self */
158e26ad7d8SSatish Balay   w3[rank]   = 0;
159e26ad7d8SSatish Balay   for (i=0; i<size; i++) {
160c1dc657dSBarry Smith     if (w1[2*i])  { w1[2*i+1] = 1; nrqs++;} /* there exists a message to proc i */
161e26ad7d8SSatish Balay   }
162b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr); /*(proc -array)*/
163e26ad7d8SSatish Balay   for (i=0,j=0; i<size; i++) {
164c1dc657dSBarry Smith     if (w1[2*i]) { pa[j] = i; j++; }
165e26ad7d8SSatish Balay   }
166e26ad7d8SSatish Balay 
167e26ad7d8SSatish Balay   /* Each message would have a header = 1 + 2*(no of IS) + data */
168e26ad7d8SSatish Balay   for (i=0; i<nrqs; i++) {
169e26ad7d8SSatish Balay     j       = pa[i];
170c1dc657dSBarry Smith     w1[2*j] += w1[2*j+1] + 2* w3[j];
171c1dc657dSBarry Smith     msz     += w1[2*j];
172e26ad7d8SSatish Balay   }
173e26ad7d8SSatish Balay   /* Do a global reduction to determine how many messages to expect*/
174c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,w1,&bsz,&nrqr);CHKERRQ(ierr);
175e26ad7d8SSatish Balay 
176e26ad7d8SSatish Balay   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
177e26ad7d8SSatish Balay   len      = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
178b0a32e0cSBarry Smith   ierr     = PetscMalloc(len,&rbuf1);CHKERRQ(ierr);
179e26ad7d8SSatish Balay   rbuf1[0] = (int*)(rbuf1 + nrqr);
180e26ad7d8SSatish Balay   for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
181e26ad7d8SSatish Balay 
182e26ad7d8SSatish Balay   /* Post the receives */
183b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&r_waits1);CHKERRQ(ierr);
184e26ad7d8SSatish Balay   for (i=0; i<nrqr; ++i) {
185e26ad7d8SSatish Balay     ierr = MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr);
186e26ad7d8SSatish Balay   }
187e26ad7d8SSatish Balay 
188e26ad7d8SSatish Balay   /* Allocate Memory for outgoing messages */
189e26ad7d8SSatish Balay   len   = 2*size*sizeof(int*) + 2*msz*sizeof(int)+ size*sizeof(int);
190b0a32e0cSBarry Smith   ierr  = PetscMalloc(len,&sbuf1);CHKERRQ(ierr);
191e26ad7d8SSatish Balay   ptr   = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
192549d3d68SSatish Balay   ierr  = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr);
193e26ad7d8SSatish Balay   /* allocate memory for outgoing data + buf to receive the first reply */
194e26ad7d8SSatish Balay   tmp   = (int*)(ptr + size);
195e26ad7d8SSatish Balay   ctr   = tmp + 2*msz;
196e26ad7d8SSatish Balay 
197e26ad7d8SSatish Balay   {
198e26ad7d8SSatish Balay     int *iptr = tmp,ict = 0;
199e26ad7d8SSatish Balay     for (i=0; i<nrqs; i++) {
200e26ad7d8SSatish Balay       j         = pa[i];
201e26ad7d8SSatish Balay       iptr     += ict;
202e26ad7d8SSatish Balay       sbuf1[j]  = iptr;
203c1dc657dSBarry Smith       ict       = w1[2*j];
204e26ad7d8SSatish Balay     }
205e26ad7d8SSatish Balay   }
206e26ad7d8SSatish Balay 
207e26ad7d8SSatish Balay   /* Form the outgoing messages */
208e26ad7d8SSatish Balay   /* Initialize the header space */
209e26ad7d8SSatish Balay   for (i=0; i<nrqs; i++) {
210e26ad7d8SSatish Balay     j           = pa[i];
211e26ad7d8SSatish Balay     sbuf1[j][0] = 0;
212549d3d68SSatish Balay     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr);
213e26ad7d8SSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
214e26ad7d8SSatish Balay   }
215e26ad7d8SSatish Balay 
216e26ad7d8SSatish Balay   /* Parse the isrow and copy data into outbuf */
217e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
218549d3d68SSatish Balay     ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr);
219e26ad7d8SSatish Balay     irow_i = irow[i];
220e26ad7d8SSatish Balay     jmax   = nrow[i];
221e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
222e26ad7d8SSatish Balay       row  = irow_i[j];
223e26ad7d8SSatish Balay       proc = rtable[row];
224e26ad7d8SSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
225e26ad7d8SSatish Balay         ctr[proc]++;
226e26ad7d8SSatish Balay         *ptr[proc] = row;
227e26ad7d8SSatish Balay         ptr[proc]++;
228e26ad7d8SSatish Balay       }
229e26ad7d8SSatish Balay     }
230e26ad7d8SSatish Balay     /* Update the headers for the current IS */
231e26ad7d8SSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
23206ef35abSSatish Balay       if ((ctr_j = ctr[j])) {
233e26ad7d8SSatish Balay         sbuf1_j        = sbuf1[j];
234e26ad7d8SSatish Balay         k              = ++sbuf1_j[0];
235e26ad7d8SSatish Balay         sbuf1_j[2*k]   = ctr_j;
236e26ad7d8SSatish Balay         sbuf1_j[2*k-1] = i;
237e26ad7d8SSatish Balay       }
238e26ad7d8SSatish Balay     }
239e26ad7d8SSatish Balay   }
240e26ad7d8SSatish Balay 
241e26ad7d8SSatish Balay   /*  Now  post the sends */
242b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
243e26ad7d8SSatish Balay   for (i=0; i<nrqs; ++i) {
244e26ad7d8SSatish Balay     j = pa[i];
245c1dc657dSBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[2*j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
246e26ad7d8SSatish Balay   }
247e26ad7d8SSatish Balay 
2484b3ec233SSatish Balay   /* Post recieves to capture the row_data from other procs */
249b0a32e0cSBarry Smith   ierr  = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr);
25087828ca2SBarry Smith   ierr  = PetscMalloc((nrqs+1)*sizeof(PetscScalar*),&rbuf2);CHKERRQ(ierr);
2514b3ec233SSatish Balay   for (i=0; i<nrqs; i++) {
252e26ad7d8SSatish Balay     j        = pa[i];
253c1dc657dSBarry Smith     count    = (w1[2*j] - (2*sbuf1[j][0] + 1))*N;
25487828ca2SBarry Smith     ierr     = PetscMalloc((count+1)*sizeof(PetscScalar),&rbuf2[i]);CHKERRQ(ierr);
2554b3ec233SSatish Balay     ierr     = MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
256e26ad7d8SSatish Balay   }
257e26ad7d8SSatish Balay 
2584b3ec233SSatish Balay   /* Receive messages(row_nos) and then, pack and send off the rowvalues
2594b3ec233SSatish Balay      to the correct processors */
260e26ad7d8SSatish Balay 
261b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
262b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr);
26387828ca2SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscScalar*),&sbuf2);CHKERRQ(ierr);
264e26ad7d8SSatish Balay 
265e26ad7d8SSatish Balay   {
26687828ca2SBarry Smith     PetscScalar *sbuf2_i,*v_start;
2674b3ec233SSatish Balay     int         s_proc;
268e26ad7d8SSatish Balay     for (i=0; i<nrqr; ++i) {
269999d9058SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
2704b3ec233SSatish Balay       s_proc          = r_status1[i].MPI_SOURCE; /* send processor */
271999d9058SBarry Smith       rbuf1_i         = rbuf1[idex]; /* Actual message from s_proc */
272999d9058SBarry Smith       /* no of rows = end - start; since start is array idex[], 0idex, whel end
273999d9058SBarry Smith          is length of the buffer - which is 1idex */
274e26ad7d8SSatish Balay       start           = 2*rbuf1_i[0] + 1;
275b0a32e0cSBarry Smith       ierr            = MPI_Get_count(r_status1+i,MPI_INT,&end);CHKERRQ(ierr);
2764b3ec233SSatish Balay       /* allocate memory sufficinet to hold all the row values */
277999d9058SBarry Smith       ierr = PetscMalloc((end-start)*N*sizeof(PetscScalar),&sbuf2[idex]);CHKERRQ(ierr);
278999d9058SBarry Smith       sbuf2_i      = sbuf2[idex];
2794b3ec233SSatish Balay       /* Now pack the data */
280e26ad7d8SSatish Balay       for (j=start; j<end; j++) {
2814b3ec233SSatish Balay         row = rbuf1_i[j] - rstart;
2824b3ec233SSatish Balay         v_start = a->v + row;
2834b3ec233SSatish Balay         for (k=0; k<N; k++) {
2844b3ec233SSatish Balay           sbuf2_i[0] = v_start[0];
285273d9f13SBarry Smith           sbuf2_i++; v_start += C->m;
286e26ad7d8SSatish Balay         }
287e26ad7d8SSatish Balay       }
2884b3ec233SSatish Balay       /* Now send off the data */
289999d9058SBarry Smith       ierr = MPI_Isend(sbuf2[idex],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);CHKERRQ(ierr);
2904b3ec233SSatish Balay     }
2914b3ec233SSatish Balay   }
2924b3ec233SSatish Balay   /* End Send-Recv of IS + row_numbers */
293606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
294606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
295b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr);
296e26ad7d8SSatish Balay   ierr      = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);
297606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
298606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
299e26ad7d8SSatish Balay 
3004b3ec233SSatish Balay   /* Create the submatrices */
3014b3ec233SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
302e26ad7d8SSatish Balay     for (i=0; i<ismax; i++) {
3034b3ec233SSatish Balay       mat = (Mat_SeqDense *)(submats[i]->data);
304273d9f13SBarry Smith       if ((submats[i]->m != nrow[i]) || (submats[i]->n != ncol[i])) {
30529bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
306e26ad7d8SSatish Balay       }
30787828ca2SBarry Smith       ierr = PetscMemzero(mat->v,submats[i]->m*submats[i]->n*sizeof(PetscScalar));CHKERRQ(ierr);
3084b3ec233SSatish Balay       submats[i]->factor = C->factor;
309e26ad7d8SSatish Balay     }
3104b3ec233SSatish Balay   } else {
311e26ad7d8SSatish Balay     for (i=0; i<ismax; i++) {
3125c5985e7SKris Buschelman       ierr = MatCreate(PETSC_COMM_SELF,nrow[i],ncol[i],nrow[i],ncol[i],submats+i);CHKERRQ(ierr);
3135c5985e7SKris Buschelman       ierr = MatSetType(submats[i],A->type_name);CHKERRQ(ierr);
3145c5985e7SKris Buschelman       ierr = MatSeqDenseSetPreallocation(submats[i],PETSC_NULL);CHKERRQ(ierr);
3154b3ec233SSatish Balay     }
3164b3ec233SSatish Balay   }
3174b3ec233SSatish Balay 
3184b3ec233SSatish Balay   /* Assemble the matrices */
3194b3ec233SSatish Balay   {
3204b3ec233SSatish Balay     int         col;
32187828ca2SBarry Smith     PetscScalar *imat_v,*mat_v,*imat_vi,*mat_vi;
3224b3ec233SSatish Balay 
3234b3ec233SSatish Balay     for (i=0; i<ismax; i++) {
3244b3ec233SSatish Balay       mat       = (Mat_SeqDense*)submats[i]->data;
3254b3ec233SSatish Balay       mat_v     = a->v;
3264b3ec233SSatish Balay       imat_v    = mat->v;
327e26ad7d8SSatish Balay       irow_i    = irow[i];
3284b3ec233SSatish Balay       m         = nrow[i];
3294b3ec233SSatish Balay       for (j=0; j<m; j++) {
330e26ad7d8SSatish Balay         row      = irow_i[j] ;
331e26ad7d8SSatish Balay         proc     = rtable[row];
332e26ad7d8SSatish Balay         if (proc == rank) {
3334b3ec233SSatish Balay           row      = row - rstart;
3344b3ec233SSatish Balay           mat_vi   = mat_v + row;
3354b3ec233SSatish Balay           imat_vi  = imat_v + j;
3364b3ec233SSatish Balay           for (k=0; k<ncol[i]; k++) {
3374b3ec233SSatish Balay             col = icol[i][k];
338273d9f13SBarry Smith             imat_vi[k*m] = mat_vi[col*C->m];
339e26ad7d8SSatish Balay           }
3404b3ec233SSatish Balay         }
341e26ad7d8SSatish Balay       }
342e26ad7d8SSatish Balay     }
343e26ad7d8SSatish Balay   }
344e26ad7d8SSatish Balay 
3454b3ec233SSatish Balay   /* Create row map. This maps c->row to submat->row for each submat*/
3464b3ec233SSatish Balay   /* this is a very expensive operation wrt memory usage */
347273d9f13SBarry Smith   len     = (1+ismax)*sizeof(int*)+ ismax*C->M*sizeof(int);
348b0a32e0cSBarry Smith   ierr    = PetscMalloc(len,&rmap);CHKERRQ(ierr);
349e26ad7d8SSatish Balay   rmap[0] = (int*)(rmap + ismax);
350273d9f13SBarry Smith   ierr    = PetscMemzero(rmap[0],ismax*C->M*sizeof(int));CHKERRQ(ierr);
351273d9f13SBarry Smith   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->M;}
352e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
353e26ad7d8SSatish Balay     rmap_i = rmap[i];
354e26ad7d8SSatish Balay     irow_i = irow[i];
355e26ad7d8SSatish Balay     jmax   = nrow[i];
356e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {
357e26ad7d8SSatish Balay       rmap_i[irow_i[j]] = j;
358e26ad7d8SSatish Balay     }
359e26ad7d8SSatish Balay   }
360e26ad7d8SSatish Balay 
3614b3ec233SSatish Balay   /* Now Receive the row_values and assemble the rest of the matrix */
362b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr);
3634b3ec233SSatish Balay 
3644b3ec233SSatish Balay   {
3654b3ec233SSatish Balay     int         is_max,tmp1,col,*sbuf1_i,is_sz;
36687828ca2SBarry Smith     PetscScalar *rbuf2_i,*imat_v,*imat_vi;
3674b3ec233SSatish Balay 
3684b3ec233SSatish Balay     for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */
3694b3ec233SSatish Balay       ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr);
3704b3ec233SSatish Balay       /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
3714b3ec233SSatish Balay       sbuf1_i = sbuf1[pa[i]];
3724b3ec233SSatish Balay       is_max  = sbuf1_i[0];
3734b3ec233SSatish Balay       ct1     = 2*is_max+1;
374e26ad7d8SSatish Balay       rbuf2_i = rbuf2[i];
3754b3ec233SSatish Balay       for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */
376e26ad7d8SSatish Balay         is_no     = sbuf1_i[2*j-1];
3774b3ec233SSatish Balay         is_sz     = sbuf1_i[2*j];
378e26ad7d8SSatish Balay         mat       = (Mat_SeqDense*)submats[is_no]->data;
3794b3ec233SSatish Balay         imat_v    = mat->v;
3804b3ec233SSatish Balay         rmap_i    = rmap[is_no];
3814b3ec233SSatish Balay         m         = nrow[is_no];
3824b3ec233SSatish Balay         for (k=0; k<is_sz; k++,rbuf2_i+=N) {  /* For each row */
3834b3ec233SSatish Balay           row      = sbuf1_i[ct1]; ct1++;
384e26ad7d8SSatish Balay           row      = rmap_i[row];
3854b3ec233SSatish Balay           imat_vi  = imat_v + row;
3864b3ec233SSatish Balay           for (l=0; l<ncol[is_no]; l++) { /* For each col */
3874b3ec233SSatish Balay             col = icol[is_no][l];
3884b3ec233SSatish Balay             imat_vi[l*m] = rbuf2_i[col];
389e26ad7d8SSatish Balay           }
390e26ad7d8SSatish Balay         }
391e26ad7d8SSatish Balay       }
392e26ad7d8SSatish Balay     }
3934b3ec233SSatish Balay   }
3944b3ec233SSatish Balay   /* End Send-Recv of row_values */
395606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
396606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
397b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr);
3984b3ec233SSatish Balay   ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);
399606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
400606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
401e26ad7d8SSatish Balay 
402e26ad7d8SSatish Balay   /* Restore the indices */
403e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
404e26ad7d8SSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
405e26ad7d8SSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
406e26ad7d8SSatish Balay   }
407e26ad7d8SSatish Balay 
408e26ad7d8SSatish Balay   /* Destroy allocated memory */
409606d414cSSatish Balay   ierr = PetscFree(irow);CHKERRQ(ierr);
410606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
411606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
412e26ad7d8SSatish Balay 
4134b3ec233SSatish Balay 
4144b3ec233SSatish Balay   for (i=0; i<nrqs; ++i) {
415606d414cSSatish Balay     ierr = PetscFree(rbuf2[i]);CHKERRQ(ierr);
4164b3ec233SSatish Balay   }
417606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
418606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
419606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
4204b3ec233SSatish Balay 
421e26ad7d8SSatish Balay   for (i=0; i<nrqr; ++i) {
422606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
423e26ad7d8SSatish Balay   }
424e26ad7d8SSatish Balay 
425606d414cSSatish Balay   ierr = PetscFree(sbuf2);CHKERRQ(ierr);
426606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
427e26ad7d8SSatish Balay 
428e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
429e26ad7d8SSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
430e26ad7d8SSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
431e26ad7d8SSatish Balay   }
432e26ad7d8SSatish Balay 
433e26ad7d8SSatish Balay   PetscFunctionReturn(0);
434e26ad7d8SSatish Balay }
435e26ad7d8SSatish Balay 
436