xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
2*be1d678aSKris Buschelman 
30500aeceSBarry Smith /*
40500aeceSBarry Smith    Support for the parallel dense matrix vector multiply
50500aeceSBarry Smith */
670f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
70500aeceSBarry Smith 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIDense"
10dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIDense(Mat mat)
110500aeceSBarry Smith {
120500aeceSBarry Smith   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
13dfbe8321SBarry Smith   PetscErrorCode ierr;
14b9b97703SBarry Smith   IS           from,to;
150500aeceSBarry Smith   Vec          gvec;
160500aeceSBarry Smith 
173a40ed3dSBarry Smith   PetscFunctionBegin;
180500aeceSBarry Smith   /* Create local vector that is used to scatter into */
19273d9f13SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,mat->N,&mdn->lvec);CHKERRQ(ierr);
200500aeceSBarry Smith 
210500aeceSBarry Smith   /* Create temporary index set for building scatter gather */
22273d9f13SBarry Smith   ierr = ISCreateStride(mat->comm,mat->N,0,1,&from);CHKERRQ(ierr);
23273d9f13SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,mat->N,0,1,&to);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 
28273d9f13SBarry Smith   ierr = VecCreateMPI(mat->comm,mdn->nvec,mat->N,&gvec);CHKERRQ(ierr);
29e19c7942SLois Curfman McInnes 
300500aeceSBarry Smith   /* Generate the scatter context */
31b9b97703SBarry Smith   ierr = VecScatterCreate(gvec,from,mdn->lvec,to,&mdn->Mvctx);CHKERRQ(ierr);
3252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,mdn->Mvctx);CHKERRQ(ierr);
3352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,mdn->lvec);CHKERRQ(ierr);
3452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
3552e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
3652e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,gvec);CHKERRQ(ierr);
37e19c7942SLois Curfman McInnes 
38b9b97703SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
39b9b97703SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
40d9f96c7cSLois Curfman McInnes   ierr = VecDestroy(gvec);CHKERRQ(ierr);
413a40ed3dSBarry Smith   PetscFunctionReturn(0);
420500aeceSBarry Smith }
430500aeceSBarry Smith 
4413f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_MPIDense_Local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*);
454a2ae208SSatish Balay #undef __FUNCT__
464a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIDense"
4713f74950SBarry Smith PetscErrorCode MatGetSubMatrices_MPIDense(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
48e26ad7d8SSatish Balay {
496849ba73SBarry Smith   PetscErrorCode ierr;
5013f74950SBarry Smith   PetscInt           nmax,nstages_local,nstages,i,pos,max_no;
51e26ad7d8SSatish Balay 
52e26ad7d8SSatish Balay   PetscFunctionBegin;
53e26ad7d8SSatish Balay   /* Allocate memory to hold all the submatrices */
54e26ad7d8SSatish Balay   if (scall != MAT_REUSE_MATRIX) {
55b0a32e0cSBarry Smith     ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr);
56e26ad7d8SSatish Balay   }
57e26ad7d8SSatish Balay   /* Determine the number of stages through which submatrices are done */
5813f74950SBarry Smith   nmax          = 20*1000000 / (C->N * sizeof(PetscInt));
59329f5518SBarry Smith   if (!nmax) nmax = 1;
60e26ad7d8SSatish Balay   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
61e26ad7d8SSatish Balay 
62e26ad7d8SSatish Balay   /* Make sure every processor loops through the nstages */
6313f74950SBarry Smith   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,C->comm);CHKERRQ(ierr);
640500aeceSBarry Smith 
650500aeceSBarry Smith 
66e26ad7d8SSatish Balay   for (i=0,pos=0; i<nstages; i++) {
67e26ad7d8SSatish Balay     if (pos+nmax <= ismax) max_no = nmax;
68e26ad7d8SSatish Balay     else if (pos == ismax) max_no = 0;
69e26ad7d8SSatish Balay     else                   max_no = ismax-pos;
70e26ad7d8SSatish Balay     ierr = MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);CHKERRQ(ierr);
71e26ad7d8SSatish Balay     pos += max_no;
72e26ad7d8SSatish Balay   }
73e26ad7d8SSatish Balay   PetscFunctionReturn(0);
74e26ad7d8SSatish Balay }
75e26ad7d8SSatish Balay /* -------------------------------------------------------------------------*/
764a2ae208SSatish Balay #undef __FUNCT__
774a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIDense_Local"
7813f74950SBarry Smith PetscErrorCode MatGetSubMatrices_MPIDense_Local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
79e26ad7d8SSatish Balay {
80e26ad7d8SSatish Balay   Mat_MPIDense   *c = (Mat_MPIDense*)C->data;
81e26ad7d8SSatish Balay   Mat            A = c->A;
82e26ad7d8SSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*mat;
836849ba73SBarry Smith   PetscErrorCode ierr;
8413f74950SBarry Smith   PetscMPIInt    rank,size,tag0,tag1,idex,end,i;
8513f74950SBarry Smith   PetscInt       N = C->N,rstart = c->rstart,count;
8613f74950SBarry Smith   PetscInt       **irow,**icol,*nrow,*ncol,*w1,*w3,*w4,*rtable,start;
8713f74950SBarry Smith   PetscInt       **sbuf1,m,j,k,l,ct1,**rbuf1,row,proc;
8813f74950SBarry Smith   PetscInt       nrqs,msz,**ptr,*ctr,*pa,*tmp,bsz,nrqr;
8913f74950SBarry Smith   PetscInt       is_no,jmax,*irow_i,**rmap,*rmap_i;
9013f74950SBarry Smith   PetscInt       len,ctr_j,*sbuf1_j,*rbuf1_i;
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;
9487828ca2SBarry Smith   PetscScalar    **rbuf2,**sbuf2;
955c5985e7SKris Buschelman   PetscTruth     sorted;
96e26ad7d8SSatish Balay 
97e26ad7d8SSatish Balay   PetscFunctionBegin;
98e26ad7d8SSatish Balay   comm   = C->comm;
99e26ad7d8SSatish Balay   tag0   = C->tag;
100e26ad7d8SSatish Balay   size   = c->size;
101e26ad7d8SSatish Balay   rank   = c->rank;
102273d9f13SBarry Smith   m      = C->M;
103e26ad7d8SSatish Balay 
104e26ad7d8SSatish Balay   /* Get some new tags to keep the communication clean */
105e26ad7d8SSatish Balay   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr);
106e26ad7d8SSatish Balay 
107e26ad7d8SSatish Balay     /* Check if the col indices are sorted */
108e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
1095c5985e7SKris Buschelman     ierr = ISSorted(isrow[i],&sorted);CHKERRQ(ierr);
1105c5985e7SKris Buschelman     if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1115c5985e7SKris Buschelman     ierr = ISSorted(iscol[i],&sorted);CHKERRQ(ierr);
1125c5985e7SKris Buschelman     if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
113e26ad7d8SSatish Balay   }
114e26ad7d8SSatish Balay 
11513f74950SBarry Smith   len    =  2*ismax*(sizeof(PetscInt*)+sizeof(PetscInt)) + (m+1)*sizeof(PetscInt);
116b0a32e0cSBarry Smith   ierr   = PetscMalloc(len,&irow);CHKERRQ(ierr);
117e26ad7d8SSatish Balay   icol   = irow + ismax;
11813f74950SBarry Smith   nrow   = (PetscInt*)(icol + ismax);
119e26ad7d8SSatish Balay   ncol   = nrow + ismax;
120e26ad7d8SSatish Balay   rtable = ncol + ismax;
121e26ad7d8SSatish Balay 
122e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
123e26ad7d8SSatish Balay     ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr);
124e26ad7d8SSatish Balay     ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr);
125b9b97703SBarry Smith     ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr);
126b9b97703SBarry Smith     ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr);
127e26ad7d8SSatish Balay   }
128e26ad7d8SSatish Balay 
129e26ad7d8SSatish Balay   /* Create hash table for the mapping :row -> proc*/
130e26ad7d8SSatish Balay   for (i=0,j=0; i<size; i++) {
131e26ad7d8SSatish Balay     jmax = c->rowners[i+1];
132e26ad7d8SSatish Balay     for (; j<jmax; j++) {
133e26ad7d8SSatish Balay       rtable[j] = i;
134e26ad7d8SSatish Balay     }
135e26ad7d8SSatish Balay   }
136e26ad7d8SSatish Balay 
137e26ad7d8SSatish Balay   /* evaluate communication - mesg to who,length of mesg, and buffer space
138e26ad7d8SSatish Balay      required. Based on this, buffers are allocated, and data copied into them*/
13913f74950SBarry Smith   ierr   = PetscMalloc(size*4*sizeof(PetscInt),&w1);CHKERRQ(ierr); /* mesg size */
140c1dc657dSBarry Smith   w3     = w1 + 2*size;      /* no of IS that needs to be sent to proc i */
141c1dc657dSBarry Smith   w4     = w3 + size;      /* temp work space used in determining w1,  w3 */
14213f74950SBarry Smith   ierr = PetscMemzero(w1,size*3*sizeof(PetscInt));CHKERRQ(ierr); /* initialize work vector*/
143e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
14413f74950SBarry Smith     ierr   = PetscMemzero(w4,size*sizeof(PetscInt));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++) {
153c1dc657dSBarry Smith       if (w4[j]) { w1[2*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) */
159c1dc657dSBarry Smith   w1[2*rank] = 0;              /* no mesg sent to self */
160e26ad7d8SSatish Balay   w3[rank]   = 0;
161e26ad7d8SSatish Balay   for (i=0; i<size; i++) {
162c1dc657dSBarry Smith     if (w1[2*i])  { w1[2*i+1] = 1; nrqs++;} /* there exists a message to proc i */
163e26ad7d8SSatish Balay   }
16413f74950SBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa);CHKERRQ(ierr); /*(proc -array)*/
165e26ad7d8SSatish Balay   for (i=0,j=0; i<size; i++) {
166c1dc657dSBarry Smith     if (w1[2*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];
172c1dc657dSBarry Smith     w1[2*j] += w1[2*j+1] + 2* w3[j];
173c1dc657dSBarry Smith     msz     += w1[2*j];
174e26ad7d8SSatish Balay   }
175e26ad7d8SSatish Balay   /* Do a global reduction to determine how many messages to expect*/
176c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,w1,&bsz,&nrqr);CHKERRQ(ierr);
177e26ad7d8SSatish Balay 
178e26ad7d8SSatish Balay   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
17913f74950SBarry Smith   len      = (nrqr+1)*sizeof(PetscInt*) + nrqr*bsz*sizeof(PetscInt);
180b0a32e0cSBarry Smith   ierr     = PetscMalloc(len,&rbuf1);CHKERRQ(ierr);
18113f74950SBarry Smith   rbuf1[0] = (PetscInt*)(rbuf1 + nrqr);
182e26ad7d8SSatish Balay   for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
183e26ad7d8SSatish Balay 
184e26ad7d8SSatish Balay   /* Post the receives */
185b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&r_waits1);CHKERRQ(ierr);
186e26ad7d8SSatish Balay   for (i=0; i<nrqr; ++i) {
18713f74950SBarry Smith     ierr = MPI_Irecv(rbuf1[i],bsz,MPIU_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr);
188e26ad7d8SSatish Balay   }
189e26ad7d8SSatish Balay 
190e26ad7d8SSatish Balay   /* Allocate Memory for outgoing messages */
19113f74950SBarry Smith   len   = 2*size*sizeof(PetscInt*) + 2*msz*sizeof(PetscInt)+ size*sizeof(PetscInt);
192b0a32e0cSBarry Smith   ierr  = PetscMalloc(len,&sbuf1);CHKERRQ(ierr);
193e26ad7d8SSatish Balay   ptr   = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
19413f74950SBarry Smith   ierr  = PetscMemzero(sbuf1,2*size*sizeof(PetscInt*));CHKERRQ(ierr);
195e26ad7d8SSatish Balay   /* allocate memory for outgoing data + buf to receive the first reply */
19613f74950SBarry Smith   tmp   = (PetscInt*)(ptr + size);
197e26ad7d8SSatish Balay   ctr   = tmp + 2*msz;
198e26ad7d8SSatish Balay 
199e26ad7d8SSatish Balay   {
20013f74950SBarry Smith     PetscInt *iptr = tmp,ict = 0;
201e26ad7d8SSatish Balay     for (i=0; i<nrqs; i++) {
202e26ad7d8SSatish Balay       j         = pa[i];
203e26ad7d8SSatish Balay       iptr     += ict;
204e26ad7d8SSatish Balay       sbuf1[j]  = iptr;
205c1dc657dSBarry Smith       ict       = w1[2*j];
206e26ad7d8SSatish Balay     }
207e26ad7d8SSatish Balay   }
208e26ad7d8SSatish Balay 
209e26ad7d8SSatish Balay   /* Form the outgoing messages */
210e26ad7d8SSatish Balay   /* Initialize the header space */
211e26ad7d8SSatish Balay   for (i=0; i<nrqs; i++) {
212e26ad7d8SSatish Balay     j           = pa[i];
213e26ad7d8SSatish Balay     sbuf1[j][0] = 0;
21413f74950SBarry Smith     ierr        = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr);
215e26ad7d8SSatish Balay     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
216e26ad7d8SSatish Balay   }
217e26ad7d8SSatish Balay 
218e26ad7d8SSatish Balay   /* Parse the isrow and copy data into outbuf */
219e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
22013f74950SBarry Smith     ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr);
221e26ad7d8SSatish Balay     irow_i = irow[i];
222e26ad7d8SSatish Balay     jmax   = nrow[i];
223e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
224e26ad7d8SSatish Balay       row  = irow_i[j];
225e26ad7d8SSatish Balay       proc = rtable[row];
226e26ad7d8SSatish Balay       if (proc != rank) { /* copy to the outgoing buf*/
227e26ad7d8SSatish Balay         ctr[proc]++;
228e26ad7d8SSatish Balay         *ptr[proc] = row;
229e26ad7d8SSatish Balay         ptr[proc]++;
230e26ad7d8SSatish Balay       }
231e26ad7d8SSatish Balay     }
232e26ad7d8SSatish Balay     /* Update the headers for the current IS */
233e26ad7d8SSatish Balay     for (j=0; j<size; j++) { /* Can Optimise this loop too */
23406ef35abSSatish Balay       if ((ctr_j = ctr[j])) {
235e26ad7d8SSatish Balay         sbuf1_j        = sbuf1[j];
236e26ad7d8SSatish Balay         k              = ++sbuf1_j[0];
237e26ad7d8SSatish Balay         sbuf1_j[2*k]   = ctr_j;
238e26ad7d8SSatish Balay         sbuf1_j[2*k-1] = i;
239e26ad7d8SSatish Balay       }
240e26ad7d8SSatish Balay     }
241e26ad7d8SSatish Balay   }
242e26ad7d8SSatish Balay 
243e26ad7d8SSatish Balay   /*  Now  post the sends */
244b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr);
245e26ad7d8SSatish Balay   for (i=0; i<nrqs; ++i) {
246e26ad7d8SSatish Balay     j = pa[i];
24713f74950SBarry Smith     ierr = MPI_Isend(sbuf1[j],w1[2*j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr);
248e26ad7d8SSatish Balay   }
249e26ad7d8SSatish Balay 
2504b3ec233SSatish Balay   /* Post recieves to capture the row_data from other procs */
251b0a32e0cSBarry Smith   ierr  = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr);
25287828ca2SBarry Smith   ierr  = PetscMalloc((nrqs+1)*sizeof(PetscScalar*),&rbuf2);CHKERRQ(ierr);
2534b3ec233SSatish Balay   for (i=0; i<nrqs; i++) {
254e26ad7d8SSatish Balay     j        = pa[i];
255c1dc657dSBarry Smith     count    = (w1[2*j] - (2*sbuf1[j][0] + 1))*N;
25687828ca2SBarry Smith     ierr     = PetscMalloc((count+1)*sizeof(PetscScalar),&rbuf2[i]);CHKERRQ(ierr);
2574b3ec233SSatish Balay     ierr     = MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);CHKERRQ(ierr);
258e26ad7d8SSatish Balay   }
259e26ad7d8SSatish Balay 
2604b3ec233SSatish Balay   /* Receive messages(row_nos) and then, pack and send off the rowvalues
2614b3ec233SSatish Balay      to the correct processors */
262e26ad7d8SSatish Balay 
263b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr);
264b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr);
26587828ca2SBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(PetscScalar*),&sbuf2);CHKERRQ(ierr);
266e26ad7d8SSatish Balay 
267e26ad7d8SSatish Balay   {
26887828ca2SBarry Smith     PetscScalar *sbuf2_i,*v_start;
26913f74950SBarry Smith     PetscInt         s_proc;
270e26ad7d8SSatish Balay     for (i=0; i<nrqr; ++i) {
271999d9058SBarry Smith       ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr);
2724b3ec233SSatish Balay       s_proc          = r_status1[i].MPI_SOURCE; /* send processor */
273999d9058SBarry Smith       rbuf1_i         = rbuf1[idex]; /* Actual message from s_proc */
274999d9058SBarry Smith       /* no of rows = end - start; since start is array idex[], 0idex, whel end
275999d9058SBarry Smith          is length of the buffer - which is 1idex */
276e26ad7d8SSatish Balay       start           = 2*rbuf1_i[0] + 1;
27713f74950SBarry Smith       ierr            = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr);
2784b3ec233SSatish Balay       /* allocate memory sufficinet to hold all the row values */
279999d9058SBarry Smith       ierr = PetscMalloc((end-start)*N*sizeof(PetscScalar),&sbuf2[idex]);CHKERRQ(ierr);
280999d9058SBarry Smith       sbuf2_i      = sbuf2[idex];
2814b3ec233SSatish Balay       /* Now pack the data */
282e26ad7d8SSatish Balay       for (j=start; j<end; j++) {
2834b3ec233SSatish Balay         row = rbuf1_i[j] - rstart;
2844b3ec233SSatish Balay         v_start = a->v + row;
2854b3ec233SSatish Balay         for (k=0; k<N; k++) {
2864b3ec233SSatish Balay           sbuf2_i[0] = v_start[0];
287273d9f13SBarry Smith           sbuf2_i++; v_start += C->m;
288e26ad7d8SSatish Balay         }
289e26ad7d8SSatish Balay       }
2904b3ec233SSatish Balay       /* Now send off the data */
291999d9058SBarry Smith       ierr = MPI_Isend(sbuf2[idex],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);CHKERRQ(ierr);
2924b3ec233SSatish Balay     }
2934b3ec233SSatish Balay   }
2944b3ec233SSatish Balay   /* End Send-Recv of IS + row_numbers */
295606d414cSSatish Balay   ierr = PetscFree(r_status1);CHKERRQ(ierr);
296606d414cSSatish Balay   ierr = PetscFree(r_waits1);CHKERRQ(ierr);
297b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr);
298e26ad7d8SSatish Balay   ierr      = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);
299606d414cSSatish Balay   ierr = PetscFree(s_status1);CHKERRQ(ierr);
300606d414cSSatish Balay   ierr = PetscFree(s_waits1);CHKERRQ(ierr);
301e26ad7d8SSatish Balay 
3024b3ec233SSatish Balay   /* Create the submatrices */
3034b3ec233SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
304e26ad7d8SSatish Balay     for (i=0; i<ismax; i++) {
3054b3ec233SSatish Balay       mat = (Mat_SeqDense *)(submats[i]->data);
306273d9f13SBarry Smith       if ((submats[i]->m != nrow[i]) || (submats[i]->n != ncol[i])) {
30729bbc08cSBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
308e26ad7d8SSatish Balay       }
30987828ca2SBarry Smith       ierr = PetscMemzero(mat->v,submats[i]->m*submats[i]->n*sizeof(PetscScalar));CHKERRQ(ierr);
3104b3ec233SSatish Balay       submats[i]->factor = C->factor;
311e26ad7d8SSatish Balay     }
3124b3ec233SSatish Balay   } else {
313e26ad7d8SSatish Balay     for (i=0; i<ismax; i++) {
3145c5985e7SKris Buschelman       ierr = MatCreate(PETSC_COMM_SELF,nrow[i],ncol[i],nrow[i],ncol[i],submats+i);CHKERRQ(ierr);
3155c5985e7SKris Buschelman       ierr = MatSetType(submats[i],A->type_name);CHKERRQ(ierr);
3165c5985e7SKris Buschelman       ierr = MatSeqDenseSetPreallocation(submats[i],PETSC_NULL);CHKERRQ(ierr);
3174b3ec233SSatish Balay     }
3184b3ec233SSatish Balay   }
3194b3ec233SSatish Balay 
3204b3ec233SSatish Balay   /* Assemble the matrices */
3214b3ec233SSatish Balay   {
32213f74950SBarry Smith     PetscInt         col;
32387828ca2SBarry Smith     PetscScalar *imat_v,*mat_v,*imat_vi,*mat_vi;
3244b3ec233SSatish Balay 
3254b3ec233SSatish Balay     for (i=0; i<ismax; i++) {
3264b3ec233SSatish Balay       mat       = (Mat_SeqDense*)submats[i]->data;
3274b3ec233SSatish Balay       mat_v     = a->v;
3284b3ec233SSatish Balay       imat_v    = mat->v;
329e26ad7d8SSatish Balay       irow_i    = irow[i];
3304b3ec233SSatish Balay       m         = nrow[i];
3314b3ec233SSatish Balay       for (j=0; j<m; j++) {
332e26ad7d8SSatish Balay         row      = irow_i[j] ;
333e26ad7d8SSatish Balay         proc     = rtable[row];
334e26ad7d8SSatish Balay         if (proc == rank) {
3354b3ec233SSatish Balay           row      = row - rstart;
3364b3ec233SSatish Balay           mat_vi   = mat_v + row;
3374b3ec233SSatish Balay           imat_vi  = imat_v + j;
3384b3ec233SSatish Balay           for (k=0; k<ncol[i]; k++) {
3394b3ec233SSatish Balay             col = icol[i][k];
340273d9f13SBarry Smith             imat_vi[k*m] = mat_vi[col*C->m];
341e26ad7d8SSatish Balay           }
3424b3ec233SSatish Balay         }
343e26ad7d8SSatish Balay       }
344e26ad7d8SSatish Balay     }
345e26ad7d8SSatish Balay   }
346e26ad7d8SSatish Balay 
3474b3ec233SSatish Balay   /* Create row map. This maps c->row to submat->row for each submat*/
3484b3ec233SSatish Balay   /* this is a very expensive operation wrt memory usage */
34913f74950SBarry Smith   len     = (1+ismax)*sizeof(PetscInt*)+ ismax*C->M*sizeof(PetscInt);
350b0a32e0cSBarry Smith   ierr    = PetscMalloc(len,&rmap);CHKERRQ(ierr);
35113f74950SBarry Smith   rmap[0] = (PetscInt*)(rmap + ismax);
35213f74950SBarry Smith   ierr    = PetscMemzero(rmap[0],ismax*C->M*sizeof(PetscInt));CHKERRQ(ierr);
353273d9f13SBarry Smith   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->M;}
354e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
355e26ad7d8SSatish Balay     rmap_i = rmap[i];
356e26ad7d8SSatish Balay     irow_i = irow[i];
357e26ad7d8SSatish Balay     jmax   = nrow[i];
358e26ad7d8SSatish Balay     for (j=0; j<jmax; j++) {
359e26ad7d8SSatish Balay       rmap_i[irow_i[j]] = j;
360e26ad7d8SSatish Balay     }
361e26ad7d8SSatish Balay   }
362e26ad7d8SSatish Balay 
3634b3ec233SSatish Balay   /* Now Receive the row_values and assemble the rest of the matrix */
364b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr);
3654b3ec233SSatish Balay 
3664b3ec233SSatish Balay   {
36713f74950SBarry Smith     PetscInt    is_max,tmp1,col,*sbuf1_i,is_sz;
36887828ca2SBarry Smith     PetscScalar *rbuf2_i,*imat_v,*imat_vi;
3694b3ec233SSatish Balay 
3704b3ec233SSatish Balay     for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */
3714b3ec233SSatish Balay       ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr);
3724b3ec233SSatish Balay       /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
3734b3ec233SSatish Balay       sbuf1_i = sbuf1[pa[i]];
3744b3ec233SSatish Balay       is_max  = sbuf1_i[0];
3754b3ec233SSatish Balay       ct1     = 2*is_max+1;
376e26ad7d8SSatish Balay       rbuf2_i = rbuf2[i];
3774b3ec233SSatish Balay       for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */
378e26ad7d8SSatish Balay         is_no     = sbuf1_i[2*j-1];
3794b3ec233SSatish Balay         is_sz     = sbuf1_i[2*j];
380e26ad7d8SSatish Balay         mat       = (Mat_SeqDense*)submats[is_no]->data;
3814b3ec233SSatish Balay         imat_v    = mat->v;
3824b3ec233SSatish Balay         rmap_i    = rmap[is_no];
3834b3ec233SSatish Balay         m         = nrow[is_no];
3844b3ec233SSatish Balay         for (k=0; k<is_sz; k++,rbuf2_i+=N) {  /* For each row */
3854b3ec233SSatish Balay           row      = sbuf1_i[ct1]; ct1++;
386e26ad7d8SSatish Balay           row      = rmap_i[row];
3874b3ec233SSatish Balay           imat_vi  = imat_v + row;
3884b3ec233SSatish Balay           for (l=0; l<ncol[is_no]; l++) { /* For each col */
3894b3ec233SSatish Balay             col = icol[is_no][l];
3904b3ec233SSatish Balay             imat_vi[l*m] = rbuf2_i[col];
391e26ad7d8SSatish Balay           }
392e26ad7d8SSatish Balay         }
393e26ad7d8SSatish Balay       }
394e26ad7d8SSatish Balay     }
3954b3ec233SSatish Balay   }
3964b3ec233SSatish Balay   /* End Send-Recv of row_values */
397606d414cSSatish Balay   ierr = PetscFree(r_status2);CHKERRQ(ierr);
398606d414cSSatish Balay   ierr = PetscFree(r_waits2);CHKERRQ(ierr);
399b0a32e0cSBarry Smith   ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr);
4004b3ec233SSatish Balay   ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);
401606d414cSSatish Balay   ierr = PetscFree(s_status2);CHKERRQ(ierr);
402606d414cSSatish Balay   ierr = PetscFree(s_waits2);CHKERRQ(ierr);
403e26ad7d8SSatish Balay 
404e26ad7d8SSatish Balay   /* Restore the indices */
405e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
406e26ad7d8SSatish Balay     ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr);
407e26ad7d8SSatish Balay     ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr);
408e26ad7d8SSatish Balay   }
409e26ad7d8SSatish Balay 
410e26ad7d8SSatish Balay   /* Destroy allocated memory */
411606d414cSSatish Balay   ierr = PetscFree(irow);CHKERRQ(ierr);
412606d414cSSatish Balay   ierr = PetscFree(w1);CHKERRQ(ierr);
413606d414cSSatish Balay   ierr = PetscFree(pa);CHKERRQ(ierr);
414e26ad7d8SSatish Balay 
4154b3ec233SSatish Balay 
4164b3ec233SSatish Balay   for (i=0; i<nrqs; ++i) {
417606d414cSSatish Balay     ierr = PetscFree(rbuf2[i]);CHKERRQ(ierr);
4184b3ec233SSatish Balay   }
419606d414cSSatish Balay   ierr = PetscFree(rbuf2);CHKERRQ(ierr);
420606d414cSSatish Balay   ierr = PetscFree(sbuf1);CHKERRQ(ierr);
421606d414cSSatish Balay   ierr = PetscFree(rbuf1);CHKERRQ(ierr);
4224b3ec233SSatish Balay 
423e26ad7d8SSatish Balay   for (i=0; i<nrqr; ++i) {
424606d414cSSatish Balay     ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr);
425e26ad7d8SSatish Balay   }
426e26ad7d8SSatish Balay 
427606d414cSSatish Balay   ierr = PetscFree(sbuf2);CHKERRQ(ierr);
428606d414cSSatish Balay   ierr = PetscFree(rmap);CHKERRQ(ierr);
429e26ad7d8SSatish Balay 
430e26ad7d8SSatish Balay   for (i=0; i<ismax; i++) {
431e26ad7d8SSatish Balay     ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
432e26ad7d8SSatish Balay     ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
433e26ad7d8SSatish Balay   }
434e26ad7d8SSatish Balay 
435e26ad7d8SSatish Balay   PetscFunctionReturn(0);
436e26ad7d8SSatish Balay }
437e26ad7d8SSatish Balay 
438