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