Actual source code: mmdense.c
1: /*$Id: mmdense.c,v 1.40 2001/09/07 20:09:22 bsmith Exp $*/
3: /*
4: Support for the parallel dense matrix vector multiply
5: */
6: #include src/mat/impls/dense/mpi/mpidense.h
10: int MatSetUpMultiply_MPIDense(Mat mat)
11: {
12: Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
13: int ierr;
14: IS from,to;
15: Vec gvec;
18: /* Create local vector that is used to scatter into */
19: VecCreateSeq(PETSC_COMM_SELF,mat->N,&mdn->lvec);
21: /* Create temporary index set for building scatter gather */
22: ISCreateStride(mat->comm,mat->N,0,1,&from);
23: ISCreateStride(PETSC_COMM_SELF,mat->N,0,1,&to);
25: /* Create temporary global vector to generate scatter context */
26: /* n = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank]; */
28: VecCreateMPI(mat->comm,mdn->nvec,mat->N,&gvec);
30: /* Generate the scatter context */
31: VecScatterCreate(gvec,from,mdn->lvec,to,&mdn->Mvctx);
32: PetscLogObjectParent(mat,mdn->Mvctx);
33: PetscLogObjectParent(mat,mdn->lvec);
34: PetscLogObjectParent(mat,from);
35: PetscLogObjectParent(mat,to);
36: PetscLogObjectParent(mat,gvec);
38: ISDestroy(to);
39: ISDestroy(from);
40: VecDestroy(gvec);
41: return(0);
42: }
44: EXTERN int MatGetSubMatrices_MPIDense_Local(Mat,int,const IS[],const IS[],MatReuse,Mat*);
47: int MatGetSubMatrices_MPIDense(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
48: {
49: int nmax,nstages_local,nstages,i,pos,max_no,ierr;
52: /* Allocate memory to hold all the submatrices */
53: if (scall != MAT_REUSE_MATRIX) {
54: PetscMalloc((ismax+1)*sizeof(Mat),submat);
55: }
56: /* Determine the number of stages through which submatrices are done */
57: nmax = 20*1000000 / (C->N * sizeof(int));
58: if (!nmax) nmax = 1;
59: nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
61: /* Make sure every processor loops through the nstages */
62: MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);
65: for (i=0,pos=0; i<nstages; i++) {
66: if (pos+nmax <= ismax) max_no = nmax;
67: else if (pos == ismax) max_no = 0;
68: else max_no = ismax-pos;
69: MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);
70: pos += max_no;
71: }
72: return(0);
73: }
74: /* -------------------------------------------------------------------------*/
77: int MatGetSubMatrices_MPIDense_Local(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
78: {
79: Mat_MPIDense *c = (Mat_MPIDense*)C->data;
80: Mat A = c->A;
81: Mat_SeqDense *a = (Mat_SeqDense*)A->data,*mat;
82: int N = C->N,rstart = c->rstart,count;
83: int **irow,**icol,*nrow,*ncol,*w1,*w3,*w4,*rtable,start,end,size;
84: int **sbuf1,rank,m,i,j,k,l,ct1,ierr,**rbuf1,row,proc;
85: int nrqs,msz,**ptr,idex,*ctr,*pa,*tmp,bsz,nrqr;
86: int is_no,jmax,*irow_i,**rmap,*rmap_i;
87: int len,ctr_j,*sbuf1_j,*rbuf1_i;
88: int tag0,tag1;
89: MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2;
90: MPI_Status *r_status1,*r_status2,*s_status1,*s_status2;
91: MPI_Comm comm;
92: PetscScalar **rbuf2,**sbuf2;
93: PetscTruth sorted;
96: comm = C->comm;
97: tag0 = C->tag;
98: size = c->size;
99: rank = c->rank;
100: m = C->M;
101:
102: /* Get some new tags to keep the communication clean */
103: PetscObjectGetNewTag((PetscObject)C,&tag1);
105: /* Check if the col indices are sorted */
106: for (i=0; i<ismax; i++) {
107: ISSorted(isrow[i],&sorted);
108: if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
109: ISSorted(iscol[i],&sorted);
110: if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
111: }
113: len = 2*ismax*(sizeof(int *)+sizeof(int)) + (m+1)*sizeof(int);
114: PetscMalloc(len,&irow);
115: icol = irow + ismax;
116: nrow = (int*)(icol + ismax);
117: ncol = nrow + ismax;
118: rtable = ncol + ismax;
120: for (i=0; i<ismax; i++) {
121: ISGetIndices(isrow[i],&irow[i]);
122: ISGetIndices(iscol[i],&icol[i]);
123: ISGetLocalSize(isrow[i],&nrow[i]);
124: ISGetLocalSize(iscol[i],&ncol[i]);
125: }
127: /* Create hash table for the mapping :row -> proc*/
128: for (i=0,j=0; i<size; i++) {
129: jmax = c->rowners[i+1];
130: for (; j<jmax; j++) {
131: rtable[j] = i;
132: }
133: }
135: /* evaluate communication - mesg to who,length of mesg, and buffer space
136: required. Based on this, buffers are allocated, and data copied into them*/
137: PetscMalloc(size*4*sizeof(int),&w1); /* mesg size */
138: w3 = w1 + 2*size; /* no of IS that needs to be sent to proc i */
139: w4 = w3 + size; /* temp work space used in determining w1, w3 */
140: PetscMemzero(w1,size*3*sizeof(int)); /* initialize work vector*/
141: for (i=0; i<ismax; i++) {
142: PetscMemzero(w4,size*sizeof(int)); /* initialize work vector*/
143: jmax = nrow[i];
144: irow_i = irow[i];
145: for (j=0; j<jmax; j++) {
146: row = irow_i[j];
147: proc = rtable[row];
148: w4[proc]++;
149: }
150: for (j=0; j<size; j++) {
151: if (w4[j]) { w1[2*j] += w4[j]; w3[j]++;}
152: }
153: }
154:
155: nrqs = 0; /* no of outgoing messages */
156: msz = 0; /* total mesg length (for all procs) */
157: w1[2*rank] = 0; /* no mesg sent to self */
158: w3[rank] = 0;
159: for (i=0; i<size; i++) {
160: if (w1[2*i]) { w1[2*i+1] = 1; nrqs++;} /* there exists a message to proc i */
161: }
162: PetscMalloc((nrqs+1)*sizeof(int),&pa); /*(proc -array)*/
163: for (i=0,j=0; i<size; i++) {
164: if (w1[2*i]) { pa[j] = i; j++; }
165: }
167: /* Each message would have a header = 1 + 2*(no of IS) + data */
168: for (i=0; i<nrqs; i++) {
169: j = pa[i];
170: w1[2*j] += w1[2*j+1] + 2* w3[j];
171: msz += w1[2*j];
172: }
173: /* Do a global reduction to determine how many messages to expect*/
174: PetscMaxSum(comm,w1,&bsz,&nrqr);
176: /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
177: len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
178: PetscMalloc(len,&rbuf1);
179: rbuf1[0] = (int*)(rbuf1 + nrqr);
180: for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
181:
182: /* Post the receives */
183: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&r_waits1);
184: for (i=0; i<nrqr; ++i) {
185: MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);
186: }
188: /* Allocate Memory for outgoing messages */
189: len = 2*size*sizeof(int*) + 2*msz*sizeof(int)+ size*sizeof(int);
190: PetscMalloc(len,&sbuf1);
191: ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */
192: PetscMemzero(sbuf1,2*size*sizeof(int*));
193: /* allocate memory for outgoing data + buf to receive the first reply */
194: tmp = (int*)(ptr + size);
195: ctr = tmp + 2*msz;
197: {
198: int *iptr = tmp,ict = 0;
199: for (i=0; i<nrqs; i++) {
200: j = pa[i];
201: iptr += ict;
202: sbuf1[j] = iptr;
203: ict = w1[2*j];
204: }
205: }
207: /* Form the outgoing messages */
208: /* Initialize the header space */
209: for (i=0; i<nrqs; i++) {
210: j = pa[i];
211: sbuf1[j][0] = 0;
212: PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));
213: ptr[j] = sbuf1[j] + 2*w3[j] + 1;
214: }
215:
216: /* Parse the isrow and copy data into outbuf */
217: for (i=0; i<ismax; i++) {
218: PetscMemzero(ctr,size*sizeof(int));
219: irow_i = irow[i];
220: jmax = nrow[i];
221: for (j=0; j<jmax; j++) { /* parse the indices of each IS */
222: row = irow_i[j];
223: proc = rtable[row];
224: if (proc != rank) { /* copy to the outgoing buf*/
225: ctr[proc]++;
226: *ptr[proc] = row;
227: ptr[proc]++;
228: }
229: }
230: /* Update the headers for the current IS */
231: for (j=0; j<size; j++) { /* Can Optimise this loop too */
232: if ((ctr_j = ctr[j])) {
233: sbuf1_j = sbuf1[j];
234: k = ++sbuf1_j[0];
235: sbuf1_j[2*k] = ctr_j;
236: sbuf1_j[2*k-1] = i;
237: }
238: }
239: }
241: /* Now post the sends */
242: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
243: for (i=0; i<nrqs; ++i) {
244: j = pa[i];
245: MPI_Isend(sbuf1[j],w1[2*j],MPI_INT,j,tag0,comm,s_waits1+i);
246: }
248: /* Post recieves to capture the row_data from other procs */
249: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);
250: PetscMalloc((nrqs+1)*sizeof(PetscScalar*),&rbuf2);
251: for (i=0; i<nrqs; i++) {
252: j = pa[i];
253: count = (w1[2*j] - (2*sbuf1[j][0] + 1))*N;
254: PetscMalloc((count+1)*sizeof(PetscScalar),&rbuf2[i]);
255: MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);
256: }
258: /* Receive messages(row_nos) and then, pack and send off the rowvalues
259: to the correct processors */
261: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);
262: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);
263: PetscMalloc((nrqr+1)*sizeof(PetscScalar*),&sbuf2);
264:
265: {
266: PetscScalar *sbuf2_i,*v_start;
267: int s_proc;
268: for (i=0; i<nrqr; ++i) {
269: MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);
270: s_proc = r_status1[i].MPI_SOURCE; /* send processor */
271: rbuf1_i = rbuf1[idex]; /* Actual message from s_proc */
272: /* no of rows = end - start; since start is array idex[], 0idex, whel end
273: is length of the buffer - which is 1idex */
274: start = 2*rbuf1_i[0] + 1;
275: MPI_Get_count(r_status1+i,MPI_INT,&end);
276: /* allocate memory sufficinet to hold all the row values */
277: PetscMalloc((end-start)*N*sizeof(PetscScalar),&sbuf2[idex]);
278: sbuf2_i = sbuf2[idex];
279: /* Now pack the data */
280: for (j=start; j<end; j++) {
281: row = rbuf1_i[j] - rstart;
282: v_start = a->v + row;
283: for (k=0; k<N; k++) {
284: sbuf2_i[0] = v_start[0];
285: sbuf2_i++; v_start += C->m;
286: }
287: }
288: /* Now send off the data */
289: MPI_Isend(sbuf2[idex],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);
290: }
291: }
292: /* End Send-Recv of IS + row_numbers */
293: PetscFree(r_status1);
294: PetscFree(r_waits1);
295: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);
296: MPI_Waitall(nrqs,s_waits1,s_status1);
297: PetscFree(s_status1);
298: PetscFree(s_waits1);
300: /* Create the submatrices */
301: if (scall == MAT_REUSE_MATRIX) {
302: for (i=0; i<ismax; i++) {
303: mat = (Mat_SeqDense *)(submats[i]->data);
304: if ((submats[i]->m != nrow[i]) || (submats[i]->n != ncol[i])) {
305: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
306: }
307: PetscMemzero(mat->v,submats[i]->m*submats[i]->n*sizeof(PetscScalar));
308: submats[i]->factor = C->factor;
309: }
310: } else {
311: for (i=0; i<ismax; i++) {
312: MatCreate(PETSC_COMM_SELF,nrow[i],ncol[i],nrow[i],ncol[i],submats+i);
313: MatSetType(submats[i],A->type_name);
314: MatSeqDenseSetPreallocation(submats[i],PETSC_NULL);
315: }
316: }
317:
318: /* Assemble the matrices */
319: {
320: int col;
321: PetscScalar *imat_v,*mat_v,*imat_vi,*mat_vi;
322:
323: for (i=0; i<ismax; i++) {
324: mat = (Mat_SeqDense*)submats[i]->data;
325: mat_v = a->v;
326: imat_v = mat->v;
327: irow_i = irow[i];
328: m = nrow[i];
329: for (j=0; j<m; j++) {
330: row = irow_i[j] ;
331: proc = rtable[row];
332: if (proc == rank) {
333: row = row - rstart;
334: mat_vi = mat_v + row;
335: imat_vi = imat_v + j;
336: for (k=0; k<ncol[i]; k++) {
337: col = icol[i][k];
338: imat_vi[k*m] = mat_vi[col*C->m];
339: }
340: }
341: }
342: }
343: }
345: /* Create row map. This maps c->row to submat->row for each submat*/
346: /* this is a very expensive operation wrt memory usage */
347: len = (1+ismax)*sizeof(int*)+ ismax*C->M*sizeof(int);
348: PetscMalloc(len,&rmap);
349: rmap[0] = (int *)(rmap + ismax);
350: PetscMemzero(rmap[0],ismax*C->M*sizeof(int));
351: for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->M;}
352: for (i=0; i<ismax; i++) {
353: rmap_i = rmap[i];
354: irow_i = irow[i];
355: jmax = nrow[i];
356: for (j=0; j<jmax; j++) {
357: rmap_i[irow_i[j]] = j;
358: }
359: }
360:
361: /* Now Receive the row_values and assemble the rest of the matrix */
362: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);
364: {
365: int is_max,tmp1,col,*sbuf1_i,is_sz;
366: PetscScalar *rbuf2_i,*imat_v,*imat_vi;
367:
368: for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */
369: MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);
370: /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
371: sbuf1_i = sbuf1[pa[i]];
372: is_max = sbuf1_i[0];
373: ct1 = 2*is_max+1;
374: rbuf2_i = rbuf2[i];
375: for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */
376: is_no = sbuf1_i[2*j-1];
377: is_sz = sbuf1_i[2*j];
378: mat = (Mat_SeqDense*)submats[is_no]->data;
379: imat_v = mat->v;
380: rmap_i = rmap[is_no];
381: m = nrow[is_no];
382: for (k=0; k<is_sz; k++,rbuf2_i+=N) { /* For each row */
383: row = sbuf1_i[ct1]; ct1++;
384: row = rmap_i[row];
385: imat_vi = imat_v + row;
386: for (l=0; l<ncol[is_no]; l++) { /* For each col */
387: col = icol[is_no][l];
388: imat_vi[l*m] = rbuf2_i[col];
389: }
390: }
391: }
392: }
393: }
394: /* End Send-Recv of row_values */
395: PetscFree(r_status2);
396: PetscFree(r_waits2);
397: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);
398: MPI_Waitall(nrqr,s_waits2,s_status2);
399: PetscFree(s_status2);
400: PetscFree(s_waits2);
402: /* Restore the indices */
403: for (i=0; i<ismax; i++) {
404: ISRestoreIndices(isrow[i],irow+i);
405: ISRestoreIndices(iscol[i],icol+i);
406: }
408: /* Destroy allocated memory */
409: PetscFree(irow);
410: PetscFree(w1);
411: PetscFree(pa);
414: for (i=0; i<nrqs; ++i) {
415: PetscFree(rbuf2[i]);
416: }
417: PetscFree(rbuf2);
418: PetscFree(sbuf1);
419: PetscFree(rbuf1);
421: for (i=0; i<nrqr; ++i) {
422: PetscFree(sbuf2[i]);
423: }
425: PetscFree(sbuf2);
426: PetscFree(rmap);
428: for (i=0; i<ismax; i++) {
429: MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);
430: MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);
431: }
433: return(0);
434: }