Actual source code: mpirowbs.c
1: #define PETSCMAT_DLL
3: #include src/mat/impls/rowbs/mpi/mpirowbs.h
5: #define CHUNCKSIZE_LOCAL 10
9: static PetscErrorCode MatFreeRowbs_Private(Mat A,int n,int *i,PetscScalar *v)
10: {
14: if (v) {
15: #if defined(PETSC_USE_LOG)
16: int len = -n*(sizeof(int)+sizeof(PetscScalar));
17: #endif
18: PetscFree(v);
19: PetscLogObjectMemory(A,len);
20: }
21: return(0);
22: }
26: static PetscErrorCode MatMallocRowbs_Private(Mat A,int n,int **i,PetscScalar **v)
27: {
29: int len;
32: if (!n) {
33: *i = 0; *v = 0;
34: } else {
35: len = n*(sizeof(int) + sizeof(PetscScalar));
36: PetscMalloc(len,v);
37: PetscLogObjectMemory(A,len);
38: *i = (int*)(*v + n);
39: }
40: return(0);
41: }
45: PetscErrorCode MatScale_MPIRowbs(Mat inA,PetscScalar alpha)
46: {
47: Mat_MPIRowbs *a = (Mat_MPIRowbs*)inA->data;
48: BSspmat *A = a->A;
49: BSsprow *vs;
50: PetscScalar *ap;
51: int i,m = inA->rmap.n,nrow,j;
55: for (i=0; i<m; i++) {
56: vs = A->rows[i];
57: nrow = vs->length;
58: ap = vs->nz;
59: for (j=0; j<nrow; j++) {
60: ap[j] *= alpha;
61: }
62: }
63: PetscLogFlops(a->nz);
64: return(0);
65: }
67: /* ----------------------------------------------------------------- */
70: static PetscErrorCode MatCreateMPIRowbs_local(Mat A,int nz,const int nnz[])
71: {
72: Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)A->data;
74: int i,len,m = A->rmap.n,*tnnz;
75: BSspmat *bsmat;
76: BSsprow *vs;
79: PetscMalloc((m+1)*sizeof(int),&tnnz);
80: if (!nnz) {
81: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
82: if (nz <= 0) nz = 1;
83: for (i=0; i<m; i++) tnnz[i] = nz;
84: nz = nz*m;
85: } else {
86: nz = 0;
87: for (i=0; i<m; i++) {
88: if (nnz[i] <= 0) tnnz[i] = 1;
89: else tnnz[i] = nnz[i];
90: nz += tnnz[i];
91: }
92: }
94: /* Allocate BlockSolve matrix context */
95: PetscNewLog(A,BSspmat,&bsif->A);
96: bsmat = bsif->A;
97: BSset_mat_icc_storage(bsmat,PETSC_FALSE);
98: BSset_mat_symmetric(bsmat,PETSC_FALSE);
99: len = m*(sizeof(BSsprow*)+ sizeof(BSsprow)) + 1;
100: PetscMalloc(len,&bsmat->rows);
101: bsmat->num_rows = m;
102: bsmat->global_num_rows = A->rmap.N;
103: bsmat->map = bsif->bsmap;
104: vs = (BSsprow*)(bsmat->rows + m);
105: for (i=0; i<m; i++) {
106: bsmat->rows[i] = vs;
107: bsif->imax[i] = tnnz[i];
108: vs->diag_ind = -1;
109: MatMallocRowbs_Private(A,tnnz[i],&(vs->col),&(vs->nz));
110: /* put zero on diagonal */
111: /*vs->length = 1;
112: vs->col[0] = i + bsif->rstart;
113: vs->nz[0] = 0.0;*/
114: vs->length = 0;
115: vs++;
116: }
117: PetscLogObjectMemory(A,sizeof(BSspmat) + len);
118: bsif->nz = 0;
119: bsif->maxnz = nz;
120: bsif->sorted = 0;
121: bsif->roworiented = PETSC_TRUE;
122: bsif->nonew = 0;
123: bsif->bs_color_single = 0;
125: PetscFree(tnnz);
126: return(0);
127: }
131: static PetscErrorCode MatSetValues_MPIRowbs_local(Mat AA,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode addv)
132: {
133: Mat_MPIRowbs *mat = (Mat_MPIRowbs*)AA->data;
134: BSspmat *A = mat->A;
135: BSsprow *vs;
137: int *rp,k,a,b,t,ii,row,nrow,i,col,l,rmax;
138: int *imax = mat->imax,nonew = mat->nonew,sorted = mat->sorted;
139: PetscScalar *ap,value;
142: for (k=0; k<m; k++) { /* loop over added rows */
143: row = im[k];
144: if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
145: if (row >= AA->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,AA->rmap.n-1);
146: vs = A->rows[row];
147: ap = vs->nz; rp = vs->col;
148: rmax = imax[row]; nrow = vs->length;
149: a = 0;
150: for (l=0; l<n; l++) { /* loop over added columns */
151: if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative col: %d",in[l]);
152: if (in[l] >= AA->cmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],AA->cmap.N-1);
153: col = in[l]; value = *v++;
154: if (!sorted) a = 0; b = nrow;
155: while (b-a > 5) {
156: t = (b+a)/2;
157: if (rp[t] > col) b = t;
158: else a = t;
159: }
160: for (i=a; i<b; i++) {
161: if (rp[i] > col) break;
162: if (rp[i] == col) {
163: if (addv == ADD_VALUES) ap[i] += value;
164: else ap[i] = value;
165: goto noinsert;
166: }
167: }
168: if (nonew) goto noinsert;
169: if (nrow >= rmax) {
170: /* there is no extra room in row, therefore enlarge */
171: int *itemp,*iout,*iin = vs->col;
172: PetscScalar *vout,*vin = vs->nz,*vtemp;
174: /* malloc new storage space */
175: imax[row] += CHUNCKSIZE_LOCAL;
176: MatMallocRowbs_Private(AA,imax[row],&itemp,&vtemp);
177: vout = vtemp; iout = itemp;
178: for (ii=0; ii<i; ii++) {
179: vout[ii] = vin[ii];
180: iout[ii] = iin[ii];
181: }
182: vout[i] = value;
183: iout[i] = col;
184: for (ii=i+1; ii<=nrow; ii++) {
185: vout[ii] = vin[ii-1];
186: iout[ii] = iin[ii-1];
187: }
188: /* free old row storage */
189: if (rmax > 0) {
190: MatFreeRowbs_Private(AA,rmax,vs->col,vs->nz);
191: }
192: vs->col = iout; vs->nz = vout;
193: rmax = imax[row];
194: mat->maxnz += CHUNCKSIZE_LOCAL;
195: mat->reallocs++;
196: } else {
197: /* shift higher columns over to make room for newie */
198: for (ii=nrow-1; ii>=i; ii--) {
199: rp[ii+1] = rp[ii];
200: ap[ii+1] = ap[ii];
201: }
202: rp[i] = col;
203: ap[i] = value;
204: }
205: nrow++;
206: mat->nz++;
207: AA->same_nonzero = PETSC_FALSE;
208: noinsert:;
209: a = i + 1;
210: }
211: vs->length = nrow;
212: }
213: return(0);
214: }
219: static PetscErrorCode MatAssemblyBegin_MPIRowbs_local(Mat A,MatAssemblyType mode)
220: {
222: return(0);
223: }
227: static PetscErrorCode MatAssemblyEnd_MPIRowbs_local(Mat AA,MatAssemblyType mode)
228: {
229: Mat_MPIRowbs *a = (Mat_MPIRowbs*)AA->data;
230: BSspmat *A = a->A;
231: BSsprow *vs;
232: int i,j,rstart = AA->rmap.rstart;
235: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
237: /* Mark location of diagonal */
238: for (i=0; i<AA->rmap.n; i++) {
239: vs = A->rows[i];
240: for (j=0; j<vs->length; j++) {
241: if (vs->col[j] == i + rstart) {
242: vs->diag_ind = j;
243: break;
244: }
245: }
246: if (vs->diag_ind == -1) {
247: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"no diagonal entry");
248: }
249: }
250: return(0);
251: }
255: static PetscErrorCode MatZeroRows_MPIRowbs_local(Mat A,PetscInt N,const PetscInt rz[],PetscScalar diag)
256: {
257: Mat_MPIRowbs *a = (Mat_MPIRowbs*)A->data;
258: BSspmat *l = a->A;
260: int i,m = A->rmap.n - 1,col,base=A->rmap.rstart;
263: if (a->keepzeroedrows) {
264: for (i=0; i<N; i++) {
265: if (rz[i] < 0 || rz[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
266: PetscMemzero(l->rows[rz[i]]->nz,l->rows[rz[i]]->length*sizeof(PetscScalar));
267: if (diag != 0.0) {
268: col=rz[i]+base;
269: MatSetValues_MPIRowbs_local(A,1,&rz[i],1,&col,&diag,INSERT_VALUES);
270: }
271: }
272: } else {
273: if (diag != 0.0) {
274: for (i=0; i<N; i++) {
275: if (rz[i] < 0 || rz[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Out of range");
276: if (l->rows[rz[i]]->length > 0) { /* in case row was completely empty */
277: l->rows[rz[i]]->length = 1;
278: l->rows[rz[i]]->nz[0] = diag;
279: l->rows[rz[i]]->col[0] = A->rmap.rstart + rz[i];
280: } else {
281: col=rz[i]+base;
282: MatSetValues_MPIRowbs_local(A,1,&rz[i],1,&col,&diag,INSERT_VALUES);
283: }
284: }
285: } else {
286: for (i=0; i<N; i++) {
287: if (rz[i] < 0 || rz[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Out of range");
288: l->rows[rz[i]]->length = 0;
289: }
290: }
291: A->same_nonzero = PETSC_FALSE;
292: }
293: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
294: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
295: return(0);
296: }
300: static PetscErrorCode MatNorm_MPIRowbs_local(Mat A,NormType type,PetscReal *norm)
301: {
302: Mat_MPIRowbs *mat = (Mat_MPIRowbs*)A->data;
303: BSsprow *vs,**rs;
304: PetscScalar *xv;
305: PetscReal sum = 0.0;
307: int *xi,nz,i,j;
310: rs = mat->A->rows;
311: if (type == NORM_FROBENIUS) {
312: for (i=0; i<A->rmap.n; i++) {
313: vs = *rs++;
314: nz = vs->length;
315: xv = vs->nz;
316: while (nz--) {
317: #if defined(PETSC_USE_COMPLEX)
318: sum += PetscRealPart(PetscConj(*xv)*(*xv)); xv++;
319: #else
320: sum += (*xv)*(*xv); xv++;
321: #endif
322: }
323: }
324: *norm = sqrt(sum);
325: } else if (type == NORM_1) { /* max column norm */
326: PetscReal *tmp;
327: PetscMalloc(A->cmap.n*sizeof(PetscReal),&tmp);
328: PetscMemzero(tmp,A->cmap.n*sizeof(PetscReal));
329: *norm = 0.0;
330: for (i=0; i<A->rmap.n; i++) {
331: vs = *rs++;
332: nz = vs->length;
333: xi = vs->col;
334: xv = vs->nz;
335: while (nz--) {
336: tmp[*xi] += PetscAbsScalar(*xv);
337: xi++; xv++;
338: }
339: }
340: for (j=0; j<A->rmap.n; j++) {
341: if (tmp[j] > *norm) *norm = tmp[j];
342: }
343: PetscFree(tmp);
344: } else if (type == NORM_INFINITY) { /* max row norm */
345: *norm = 0.0;
346: for (i=0; i<A->rmap.n; i++) {
347: vs = *rs++;
348: nz = vs->length;
349: xv = vs->nz;
350: sum = 0.0;
351: while (nz--) {
352: sum += PetscAbsScalar(*xv); xv++;
353: }
354: if (sum > *norm) *norm = sum;
355: }
356: } else {
357: SETERRQ(PETSC_ERR_SUP,"No support for the two norm");
358: }
359: return(0);
360: }
362: /* ----------------------------------------------------------------- */
366: PetscErrorCode MatSetValues_MPIRowbs(Mat mat,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode av)
367: {
368: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
370: int i,j,row,col,rstart = mat->rmap.rstart,rend = mat->rmap.rend;
371: PetscTruth roworiented = a->roworiented;
374: /* Note: There's no need to "unscale" the matrix, since scaling is
375: confined to a->pA, and we're working with a->A here */
376: for (i=0; i<m; i++) {
377: if (im[i] < 0) continue;
378: if (im[i] >= mat->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",im[i],mat->rmap.N-1);
379: if (im[i] >= rstart && im[i] < rend) {
380: row = im[i] - rstart;
381: for (j=0; j<n; j++) {
382: if (in[j] < 0) continue;
383: if (in[j] >= mat->cmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[j],mat->cmap.N-1);
384: if (in[j] >= 0 && in[j] < mat->cmap.N){
385: col = in[j];
386: if (roworiented) {
387: MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,v+i*n+j,av);
388: } else {
389: MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,v+i+j*m,av);
390: }
391: } else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid column");}
392: }
393: } else {
394: if (!a->donotstash) {
395: if (roworiented) {
396: MatStashValuesRow_Private(&mat->stash,im[i],n,in,v+i*n);
397: } else {
398: MatStashValuesCol_Private(&mat->stash,im[i],n,in,v+i,m);
399: }
400: }
401: }
402: }
403: return(0);
404: }
408: PetscErrorCode MatAssemblyBegin_MPIRowbs(Mat mat,MatAssemblyType mode)
409: {
410: MPI_Comm comm = ((PetscObject)mat)->comm;
412: int nstash,reallocs;
413: InsertMode addv;
416: /* Note: There's no need to "unscale" the matrix, since scaling is
417: confined to a->pA, and we're working with a->A here */
419: /* make sure all processors are either in INSERTMODE or ADDMODE */
420: MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
421: if (addv == (ADD_VALUES|INSERT_VALUES)) {
422: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Some procs inserted; others added");
423: }
424: mat->insertmode = addv; /* in case this processor had no cache */
426: MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap.range);
427: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
428: PetscInfo2(mat,"Block-Stash has %d entries, uses %d mallocs.\n",nstash,reallocs);
429: return(0);
430: }
434: static PetscErrorCode MatView_MPIRowbs_ASCII(Mat mat,PetscViewer viewer)
435: {
436: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
438: int i,j;
439: PetscTruth iascii;
440: BSspmat *A = a->A;
441: BSsprow **rs = A->rows;
442: PetscViewerFormat format;
445: PetscViewerGetFormat(viewer,&format);
446: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
448: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
449: int ind_l,ind_g,clq_l,clq_g,color;
450: ind_l = BSlocal_num_inodes(a->pA);CHKERRBS(0);
451: ind_g = BSglobal_num_inodes(a->pA);CHKERRBS(0);
452: clq_l = BSlocal_num_cliques(a->pA);CHKERRBS(0);
453: clq_g = BSglobal_num_cliques(a->pA);CHKERRBS(0);
454: color = BSnum_colors(a->pA);CHKERRBS(0);
455: PetscViewerASCIIPrintf(viewer," %d global inode(s), %d global clique(s), %d color(s)\n",ind_g,clq_g,color);
456: PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d local inode(s), %d local clique(s)\n",a->rank,ind_l,clq_l);
457: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
458: for (i=0; i<A->num_rows; i++) {
459: PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+mat->rmap.rstart);
460: for (j=0; j<rs[i]->length; j++) {
461: if (rs[i]->nz[j]) {PetscViewerASCIISynchronizedPrintf(viewer," %d %g ",rs[i]->col[j],rs[i]->nz[j]);}
462: }
463: PetscViewerASCIISynchronizedPrintf(viewer,"\n");
464: }
465: } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
466: SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
467: } else {
468: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
469: for (i=0; i<A->num_rows; i++) {
470: PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+mat->rmap.rstart);
471: for (j=0; j<rs[i]->length; j++) {
472: PetscViewerASCIISynchronizedPrintf(viewer," %d %g ",rs[i]->col[j],rs[i]->nz[j]);
473: }
474: PetscViewerASCIISynchronizedPrintf(viewer,"\n");
475: }
476: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
477: }
478: PetscViewerFlush(viewer);
479: return(0);
480: }
484: static PetscErrorCode MatView_MPIRowbs_Binary(Mat mat,PetscViewer viewer)
485: {
486: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
488: PetscMPIInt rank,size;
489: PetscInt i,M,m,*sbuff,*rowlengths;
490: PetscInt *recvcts,*recvdisp,fd,*cols,maxnz,nz,j;
491: BSspmat *A = a->A;
492: BSsprow **rs = A->rows;
493: MPI_Comm comm = ((PetscObject)mat)->comm;
494: MPI_Status status;
495: PetscScalar *vals;
496: MatInfo info;
499: MPI_Comm_size(comm,&size);
500: MPI_Comm_rank(comm,&rank);
502: M = mat->rmap.N; m = mat->rmap.n;
503: /* First gather together on the first processor the lengths of
504: each row, and write them out to the file */
505: PetscMalloc(m*sizeof(int),&sbuff);
506: for (i=0; i<A->num_rows; i++) {
507: sbuff[i] = rs[i]->length;
508: }
509: MatGetInfo(mat,MAT_GLOBAL_SUM,&info);
510: if (!rank) {
511: PetscViewerBinaryGetDescriptor(viewer,&fd);
512: PetscMalloc((4+M)*sizeof(int),&rowlengths);
513: PetscMalloc(size*sizeof(int),&recvcts);
514: recvdisp = mat->rmap.range;
515: for (i=0; i<size; i++) {
516: recvcts[i] = recvdisp[i+1] - recvdisp[i];
517: }
518: /* first four elements of rowlength are the header */
519: rowlengths[0] = ((PetscObject)mat)->cookie;
520: rowlengths[1] = mat->rmap.N;
521: rowlengths[2] = mat->cmap.N;
522: rowlengths[3] = (int)info.nz_used;
523: MPI_Gatherv(sbuff,m,MPI_INT,rowlengths+4,recvcts,recvdisp,MPI_INT,0,comm);
524: PetscFree(sbuff);
525: PetscBinaryWrite(fd,rowlengths,4+M,PETSC_INT,PETSC_FALSE);
526: /* count the number of nonzeros on each processor */
527: PetscMemzero(recvcts,size*sizeof(int));
528: for (i=0; i<size; i++) {
529: for (j=recvdisp[i]; j<recvdisp[i+1]; j++) {
530: recvcts[i] += rowlengths[j+3];
531: }
532: }
533: /* allocate buffer long enough to hold largest one */
534: maxnz = 0;
535: for (i=0; i<size; i++) {
536: maxnz = PetscMax(maxnz,recvcts[i]);
537: }
538: PetscFree(rowlengths);
539: PetscFree(recvcts);
540: PetscMalloc(maxnz*sizeof(int),&cols);
542: /* binary store column indices for 0th processor */
543: nz = 0;
544: for (i=0; i<A->num_rows; i++) {
545: for (j=0; j<rs[i]->length; j++) {
546: cols[nz++] = rs[i]->col[j];
547: }
548: }
549: PetscBinaryWrite(fd,cols,nz,PETSC_INT,PETSC_FALSE);
551: /* receive and store column indices for all other processors */
552: for (i=1; i<size; i++) {
553: /* should tell processor that I am now ready and to begin the send */
554: MPI_Recv(cols,maxnz,MPI_INT,i,((PetscObject)mat)->tag,comm,&status);
555: MPI_Get_count(&status,MPI_INT,&nz);
556: PetscBinaryWrite(fd,cols,nz,PETSC_INT,PETSC_FALSE);
557: }
558: PetscFree(cols);
559: PetscMalloc(maxnz*sizeof(PetscScalar),&vals);
561: /* binary store values for 0th processor */
562: nz = 0;
563: for (i=0; i<A->num_rows; i++) {
564: for (j=0; j<rs[i]->length; j++) {
565: vals[nz++] = rs[i]->nz[j];
566: }
567: }
568: PetscBinaryWrite(fd,vals,nz,PETSC_SCALAR,PETSC_FALSE);
570: /* receive and store nonzeros for all other processors */
571: for (i=1; i<size; i++) {
572: /* should tell processor that I am now ready and to begin the send */
573: MPI_Recv(vals,maxnz,MPIU_SCALAR,i,((PetscObject)mat)->tag,comm,&status);
574: MPI_Get_count(&status,MPIU_SCALAR,&nz);
575: PetscBinaryWrite(fd,vals,nz,PETSC_SCALAR,PETSC_FALSE);
576: }
577: PetscFree(vals);
578: } else {
579: MPI_Gatherv(sbuff,m,MPI_INT,0,0,0,MPI_INT,0,comm);
580: PetscFree(sbuff);
582: /* count local nonzeros */
583: nz = 0;
584: for (i=0; i<A->num_rows; i++) {
585: for (j=0; j<rs[i]->length; j++) {
586: nz++;
587: }
588: }
589: /* copy into buffer column indices */
590: PetscMalloc(nz*sizeof(int),&cols);
591: nz = 0;
592: for (i=0; i<A->num_rows; i++) {
593: for (j=0; j<rs[i]->length; j++) {
594: cols[nz++] = rs[i]->col[j];
595: }
596: }
597: /* send */ /* should wait until processor zero tells me to go */
598: MPI_Send(cols,nz,MPI_INT,0,((PetscObject)mat)->tag,comm);
599: PetscFree(cols);
601: /* copy into buffer column values */
602: PetscMalloc(nz*sizeof(PetscScalar),&vals);
603: nz = 0;
604: for (i=0; i<A->num_rows; i++) {
605: for (j=0; j<rs[i]->length; j++) {
606: vals[nz++] = rs[i]->nz[j];
607: }
608: }
609: /* send */ /* should wait until processor zero tells me to go */
610: MPI_Send(vals,nz,MPIU_SCALAR,0,((PetscObject)mat)->tag,comm);
611: PetscFree(vals);
612: }
614: return(0);
615: }
619: PetscErrorCode MatView_MPIRowbs(Mat mat,PetscViewer viewer)
620: {
621: Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)mat->data;
623: PetscTruth iascii,isbinary;
626: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
627: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
628: if (!bsif->blocksolveassembly) {
629: MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
630: }
631: if (iascii) {
632: MatView_MPIRowbs_ASCII(mat,viewer);
633: } else if (isbinary) {
634: MatView_MPIRowbs_Binary(mat,viewer);
635: } else {
636: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPIRowbs matrices",((PetscObject)viewer)->type_name);
637: }
638: return(0);
639: }
640:
643: static PetscErrorCode MatAssemblyEnd_MPIRowbs_MakeSymmetric(Mat mat)
644: {
645: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
646: BSspmat *A = a->A;
647: BSsprow *vs;
648: int size,rank,M,rstart,tag,i,j,*rtable,*w1,*w3,*w4,len,proc,nrqs;
649: int msz,*pa,bsz,nrqr,**rbuf1,**sbuf1,**ptr,*tmp,*ctr,col,idx,row;
651: int ctr_j,*sbuf1_j,k;
652: PetscScalar val=0.0;
653: MPI_Comm comm;
654: MPI_Request *s_waits1,*r_waits1;
655: MPI_Status *s_status,*r_status;
658: comm = ((PetscObject)mat)->comm;
659: tag = ((PetscObject)mat)->tag;
660: size = a->size;
661: rank = a->rank;
662: M = mat->rmap.N;
663: rstart = mat->rmap.rstart;
665: PetscMalloc(M*sizeof(int),&rtable);
666: /* Create hash table for the mapping :row -> proc */
667: for (i=0,j=0; i<size; i++) {
668: len = mat->rmap.range[i+1];
669: for (; j<len; j++) {
670: rtable[j] = i;
671: }
672: }
674: /* Evaluate communication - mesg to whom, length of mesg, and buffer space
675: required. Based on this, buffers are allocated, and data copied into them. */
676: PetscMalloc(size*4*sizeof(int),&w1);/* mesg size */
677: w3 = w1 + 2*size; /* no of IS that needs to be sent to proc i */
678: w4 = w3 + size; /* temp work space used in determining w1, w3 */
679: PetscMemzero(w1,size*3*sizeof(int)); /* initialize work vector */
681: for (i=0; i<mat->rmap.n; i++) {
682: PetscMemzero(w4,size*sizeof(int)); /* initialize work vector */
683: vs = A->rows[i];
684: for (j=0; j<vs->length; j++) {
685: proc = rtable[vs->col[j]];
686: w4[proc]++;
687: }
688: for (j=0; j<size; j++) {
689: if (w4[j]) { w1[2*j] += w4[j]; w3[j]++;}
690: }
691: }
692:
693: nrqs = 0; /* number of outgoing messages */
694: msz = 0; /* total mesg length (for all proc */
695: w1[2*rank] = 0; /* no mesg sent to itself */
696: w3[rank] = 0;
697: for (i=0; i<size; i++) {
698: if (w1[2*i]) {w1[2*i+1] = 1; nrqs++;} /* there exists a message to proc i */
699: }
700: /* pa - is list of processors to communicate with */
701: PetscMalloc((nrqs+1)*sizeof(int),&pa);
702: for (i=0,j=0; i<size; i++) {
703: if (w1[2*i]) {pa[j] = i; j++;}
704: }
706: /* Each message would have a header = 1 + 2*(no of ROWS) + data */
707: for (i=0; i<nrqs; i++) {
708: j = pa[i];
709: w1[2*j] += w1[2*j+1] + 2*w3[j];
710: msz += w1[2*j];
711: }
712:
713: /* Do a global reduction to determine how many messages to expect */
714: PetscMaxSum(comm,w1,&bsz,&nrqr);
716: /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
717: len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
718: PetscMalloc(len,&rbuf1);
719: rbuf1[0] = (int*)(rbuf1 + nrqr);
720: for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
722: /* Post the receives */
723: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&r_waits1);
724: for (i=0; i<nrqr; ++i){
725: MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits1+i);
726: }
727:
728: /* Allocate Memory for outgoing messages */
729: len = 2*size*sizeof(int*) + (size+msz)*sizeof(int);
730: PetscMalloc(len,&sbuf1);
731: ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */
732: PetscMemzero(sbuf1,2*size*sizeof(int*));
733: tmp = (int*)(sbuf1 + 2*size);
734: ctr = tmp + msz;
736: {
737: int *iptr = tmp,ict = 0;
738: for (i=0; i<nrqs; i++) {
739: j = pa[i];
740: iptr += ict;
741: sbuf1[j] = iptr;
742: ict = w1[2*j];
743: }
744: }
746: /* Form the outgoing messages */
747: /* Clean up the header space */
748: for (i=0; i<nrqs; i++) {
749: j = pa[i];
750: sbuf1[j][0] = 0;
751: PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));
752: ptr[j] = sbuf1[j] + 2*w3[j] + 1;
753: }
755: /* Parse the matrix and copy the data into sbuf1 */
756: for (i=0; i<mat->rmap.n; i++) {
757: PetscMemzero(ctr,size*sizeof(int));
758: vs = A->rows[i];
759: for (j=0; j<vs->length; j++) {
760: col = vs->col[j];
761: proc = rtable[col];
762: if (proc != rank) { /* copy to the outgoing buffer */
763: ctr[proc]++;
764: *ptr[proc] = col;
765: ptr[proc]++;
766: } else {
767: row = col - rstart;
768: col = i + rstart;
769: MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,&val,ADD_VALUES);
770: }
771: }
772: /* Update the headers for the current row */
773: for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
774: if ((ctr_j = ctr[j])) {
775: sbuf1_j = sbuf1[j];
776: k = ++sbuf1_j[0];
777: sbuf1_j[2*k] = ctr_j;
778: sbuf1_j[2*k-1] = i + rstart;
779: }
780: }
781: }
782: /* Check Validity of the outgoing messages */
783: {
784: int sum;
785: for (i=0 ; i<nrqs ; i++) {
786: j = pa[i];
787: if (w3[j] != sbuf1[j][0]) {SETERRQ(PETSC_ERR_PLIB,"Blew it! Header[1] mismatch!\n"); }
788: }
790: for (i=0 ; i<nrqs ; i++) {
791: j = pa[i];
792: sum = 1;
793: for (k = 1; k <= w3[j]; k++) sum += sbuf1[j][2*k]+2;
794: if (sum != w1[2*j]) { SETERRQ(PETSC_ERR_PLIB,"Blew it! Header[2-n] mismatch!\n"); }
795: }
796: }
797:
798: /* Now post the sends */
799: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
800: for (i=0; i<nrqs; ++i) {
801: j = pa[i];
802: MPI_Isend(sbuf1[j],w1[2*j],MPI_INT,j,tag,comm,s_waits1+i);
803: }
804:
805: /* Receive messages*/
806: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status);
807: for (i=0; i<nrqr; ++i) {
808: MPI_Waitany(nrqr,r_waits1,&idx,r_status+i);
809: /* Process the Message */
810: {
811: int *rbuf1_i,n_row,ct1;
813: rbuf1_i = rbuf1[idx];
814: n_row = rbuf1_i[0];
815: ct1 = 2*n_row+1;
816: val = 0.0;
817: /* Optimise this later */
818: for (j=1; j<=n_row; j++) {
819: col = rbuf1_i[2*j-1];
820: for (k=0; k<rbuf1_i[2*j]; k++,ct1++) {
821: row = rbuf1_i[ct1] - rstart;
822: MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,&val,ADD_VALUES);
823: }
824: }
825: }
826: }
828: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);
829: if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status);}
831: PetscFree(rtable);
832: PetscFree(w1);
833: PetscFree(pa);
834: PetscFree(rbuf1);
835: PetscFree(sbuf1);
836: PetscFree(r_waits1);
837: PetscFree(s_waits1);
838: PetscFree(r_status);
839: PetscFree(s_status);
840: return(0);
841: }
843: /*
844: This does the BlockSolve portion of the matrix assembly.
845: It is provided in a separate routine so that users can
846: operate on the matrix (using MatScale(), MatShift() etc.) after
847: the matrix has been assembled but before BlockSolve has sucked it
848: in and devoured it.
849: */
852: PetscErrorCode MatAssemblyEnd_MPIRowbs_ForBlockSolve(Mat mat)
853: {
854: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
856: int ldim,low,high,i;
857: PetscScalar *diag;
860: if ((mat->was_assembled) && (!mat->same_nonzero)) { /* Free the old info */
861: if (a->pA) {BSfree_par_mat(a->pA);CHKERRBS(0);}
862: if (a->comm_pA) {BSfree_comm(a->comm_pA);CHKERRBS(0);}
863: }
865: if ((!mat->same_nonzero) || (!mat->was_assembled)) {
866: /* Indicates bypassing cliques in coloring */
867: if (a->bs_color_single) {
868: BSctx_set_si(a->procinfo,100);
869: }
870: /* Form permuted matrix for efficient parallel execution */
871: a->pA = BSmain_perm(a->procinfo,a->A);CHKERRBS(0);
872: /* Set up the communication */
873: a->comm_pA = BSsetup_forward(a->pA,a->procinfo);CHKERRBS(0);
874: } else {
875: /* Repermute the matrix */
876: BSmain_reperm(a->procinfo,a->A,a->pA);CHKERRBS(0);
877: }
879: /* Symmetrically scale the matrix by the diagonal */
880: BSscale_diag(a->pA,a->pA->diag,a->procinfo);CHKERRBS(0);
882: /* Store inverse of square root of permuted diagonal scaling matrix */
883: VecGetLocalSize(a->diag,&ldim);
884: VecGetOwnershipRange(a->diag,&low,&high);
885: VecGetArray(a->diag,&diag);
886: for (i=0; i<ldim; i++) {
887: if (a->pA->scale_diag[i] != 0.0) {
888: diag[i] = 1.0/sqrt(PetscAbsScalar(a->pA->scale_diag[i]));
889: } else {
890: diag[i] = 1.0;
891: }
892: }
893: VecRestoreArray(a->diag,&diag);
894: a->assembled_icc_storage = a->A->icc_storage;
895: a->blocksolveassembly = 1;
896: mat->was_assembled = PETSC_TRUE;
897: mat->same_nonzero = PETSC_TRUE;
898: PetscInfo(mat,"Completed BlockSolve95 matrix assembly\n");
899: return(0);
900: }
904: PetscErrorCode MatAssemblyEnd_MPIRowbs(Mat mat,MatAssemblyType mode)
905: {
906: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
908: int i,n,row,col,*rows,*cols,rstart,nzcount,flg,j,ncols;
909: PetscScalar *vals,val;
910: InsertMode addv = mat->insertmode;
913: while (1) {
914: MatStashScatterGetMesg_Private(&mat->stash,&n,&rows,&cols,&vals,&flg);
915: if (!flg) break;
916:
917: for (i=0; i<n;) {
918: /* Now identify the consecutive vals belonging to the same row */
919: for (j=i,rstart=rows[j]; j<n; j++) { if (rows[j] != rstart) break; }
920: if (j < n) ncols = j-i;
921: else ncols = n-i;
922: /* Now assemble all these values with a single function call */
923: MatSetValues_MPIRowbs(mat,1,rows+i,ncols,cols+i,vals+i,addv);
924: i = j;
925: }
926: }
927: MatStashScatterEnd_Private(&mat->stash);
929: rstart = mat->rmap.rstart;
930: nzcount = a->nz; /* This is the number of nonzeros entered by the user */
931: /* BlockSolve requires that the matrix is structurally symmetric */
932: if (mode == MAT_FINAL_ASSEMBLY && !mat->structurally_symmetric) {
933: MatAssemblyEnd_MPIRowbs_MakeSymmetric(mat);
934: }
935:
936: /* BlockSolve requires that all the diagonal elements are set */
937: val = 0.0;
938: for (i=0; i<mat->rmap.n; i++) {
939: row = i; col = i + rstart;
940: MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,&val,ADD_VALUES);
941: }
942:
943: MatAssemblyBegin_MPIRowbs_local(mat,mode);
944: MatAssemblyEnd_MPIRowbs_local(mat,mode);
945:
946: a->blocksolveassembly = 0;
947: PetscInfo4(mat,"Matrix size: %d X %d; storage space: %d unneeded,%d used\n",mat->rmap.n,mat->cmap.n,a->maxnz-a->nz,a->nz);
948: PetscInfo2(mat,"User entered %d nonzeros, PETSc added %d\n",nzcount,a->nz-nzcount);
949: PetscInfo1(mat,"Number of mallocs during MatSetValues is %d\n",a->reallocs);
950: return(0);
951: }
955: PetscErrorCode MatZeroEntries_MPIRowbs(Mat mat)
956: {
957: Mat_MPIRowbs *l = (Mat_MPIRowbs*)mat->data;
958: BSspmat *A = l->A;
959: BSsprow *vs;
960: int i,j;
963: for (i=0; i <mat->rmap.n; i++) {
964: vs = A->rows[i];
965: for (j=0; j< vs->length; j++) vs->nz[j] = 0.0;
966: }
967: return(0);
968: }
970: /* the code does not do the diagonal entries correctly unless the
971: matrix is square and the column and row owerships are identical.
972: This is a BUG.
973: */
977: PetscErrorCode MatZeroRows_MPIRowbs(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
978: {
979: Mat_MPIRowbs *l = (Mat_MPIRowbs*)A->data;
981: int i,*owners = A->rmap.range,size = l->size;
982: int *nprocs,j,idx,nsends;
983: int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
984: int *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,n,*source;
985: int *lens,imdex,*lrows,*values;
986: MPI_Comm comm = ((PetscObject)A)->comm;
987: MPI_Request *send_waits,*recv_waits;
988: MPI_Status recv_status,*send_status;
989: PetscTruth found;
992: /* first count number of contributors to each processor */
993: PetscMalloc(2*size*sizeof(int),&nprocs);
994: PetscMemzero(nprocs,2*size*sizeof(int));
995: PetscMalloc((N+1)*sizeof(int),&owner); /* see note*/
996: for (i=0; i<N; i++) {
997: idx = rows[i];
998: found = PETSC_FALSE;
999: for (j=0; j<size; j++) {
1000: if (idx >= owners[j] && idx < owners[j+1]) {
1001: nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
1002: }
1003: }
1004: if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
1005: }
1006: nsends = 0; for (i=0; i<size; i++) {nsends += nprocs[2*i+1];}
1008: /* inform other processors of number of messages and max length*/
1009: PetscMaxSum(comm,nprocs,&nmax,&nrecvs);
1011: /* post receives: */
1012: PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);
1013: PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);
1014: for (i=0; i<nrecvs; i++) {
1015: MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
1016: }
1018: /* do sends:
1019: 1) starts[i] gives the starting index in svalues for stuff going to
1020: the ith processor
1021: */
1022: PetscMalloc((N+1)*sizeof(int),&svalues);
1023: PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);
1024: PetscMalloc((size+1)*sizeof(int),&starts);
1025: starts[0] = 0;
1026: for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
1027: for (i=0; i<N; i++) {
1028: svalues[starts[owner[i]]++] = rows[i];
1029: }
1031: starts[0] = 0;
1032: for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
1033: count = 0;
1034: for (i=0; i<size; i++) {
1035: if (nprocs[2*i+1]) {
1036: MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);
1037: }
1038: }
1039: PetscFree(starts);
1041: base = owners[rank];
1043: /* wait on receives */
1044: PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);
1045: source = lens + nrecvs;
1046: count = nrecvs; slen = 0;
1047: while (count) {
1048: MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
1049: /* unpack receives into our local space */
1050: MPI_Get_count(&recv_status,MPI_INT,&n);
1051: source[imdex] = recv_status.MPI_SOURCE;
1052: lens[imdex] = n;
1053: slen += n;
1054: count--;
1055: }
1056: PetscFree(recv_waits);
1057:
1058: /* move the data into the send scatter */
1059: PetscMalloc((slen+1)*sizeof(int),&lrows);
1060: count = 0;
1061: for (i=0; i<nrecvs; i++) {
1062: values = rvalues + i*nmax;
1063: for (j=0; j<lens[i]; j++) {
1064: lrows[count++] = values[j] - base;
1065: }
1066: }
1067: PetscFree(rvalues);
1068: PetscFree(lens);
1069: PetscFree(owner);
1070: PetscFree(nprocs);
1071:
1072: /* actually zap the local rows */
1073: MatZeroRows_MPIRowbs_local(A,slen,lrows,diag);
1074: PetscFree(lrows);
1076: /* wait on sends */
1077: if (nsends) {
1078: PetscMalloc(nsends*sizeof(MPI_Status),&send_status);
1079: MPI_Waitall(nsends,send_waits,send_status);
1080: PetscFree(send_status);
1081: }
1082: PetscFree(send_waits);
1083: PetscFree(svalues);
1085: return(0);
1086: }
1090: PetscErrorCode MatNorm_MPIRowbs(Mat mat,NormType type,PetscReal *norm)
1091: {
1092: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
1093: BSsprow *vs,**rs;
1094: PetscScalar *xv;
1095: PetscReal sum = 0.0;
1097: int *xi,nz,i,j;
1100: if (a->size == 1) {
1101: MatNorm_MPIRowbs_local(mat,type,norm);
1102: } else {
1103: rs = a->A->rows;
1104: if (type == NORM_FROBENIUS) {
1105: for (i=0; i<mat->rmap.n; i++) {
1106: vs = *rs++;
1107: nz = vs->length;
1108: xv = vs->nz;
1109: while (nz--) {
1110: #if defined(PETSC_USE_COMPLEX)
1111: sum += PetscRealPart(PetscConj(*xv)*(*xv)); xv++;
1112: #else
1113: sum += (*xv)*(*xv); xv++;
1114: #endif
1115: }
1116: }
1117: MPI_Allreduce(&sum,norm,1,MPIU_REAL,MPI_SUM,((PetscObject)mat)->comm);
1118: *norm = sqrt(*norm);
1119: } else if (type == NORM_1) { /* max column norm */
1120: PetscReal *tmp,*tmp2;
1121: PetscMalloc(mat->cmap.n*sizeof(PetscReal),&tmp);
1122: PetscMalloc(mat->cmap.n*sizeof(PetscReal),&tmp2);
1123: PetscMemzero(tmp,mat->cmap.n*sizeof(PetscReal));
1124: *norm = 0.0;
1125: for (i=0; i<mat->rmap.n; i++) {
1126: vs = *rs++;
1127: nz = vs->length;
1128: xi = vs->col;
1129: xv = vs->nz;
1130: while (nz--) {
1131: tmp[*xi] += PetscAbsScalar(*xv);
1132: xi++; xv++;
1133: }
1134: }
1135: MPI_Allreduce(tmp,tmp2,mat->cmap.N,MPIU_REAL,MPI_SUM,((PetscObject)mat)->comm);
1136: for (j=0; j<mat->cmap.n; j++) {
1137: if (tmp2[j] > *norm) *norm = tmp2[j];
1138: }
1139: PetscFree(tmp);
1140: PetscFree(tmp2);
1141: } else if (type == NORM_INFINITY) { /* max row norm */
1142: PetscReal ntemp = 0.0;
1143: for (i=0; i<mat->rmap.n; i++) {
1144: vs = *rs++;
1145: nz = vs->length;
1146: xv = vs->nz;
1147: sum = 0.0;
1148: while (nz--) {
1149: sum += PetscAbsScalar(*xv); xv++;
1150: }
1151: if (sum > ntemp) ntemp = sum;
1152: }
1153: MPI_Allreduce(&ntemp,norm,1,MPIU_REAL,MPI_MAX,((PetscObject)mat)->comm);
1154: } else {
1155: SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1156: }
1157: }
1158: return(0);
1159: }
1163: PetscErrorCode MatMult_MPIRowbs(Mat mat,Vec xx,Vec yy)
1164: {
1165: Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)mat->data;
1166: BSprocinfo *bspinfo = bsif->procinfo;
1167: PetscScalar *xxa,*xworka,*yya;
1171: if (!bsif->blocksolveassembly) {
1172: MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
1173: }
1175: /* Permute and apply diagonal scaling: [ xwork = D^{1/2} * x ] */
1176: if (!bsif->vecs_permscale) {
1177: VecGetArray(bsif->xwork,&xworka);
1178: VecGetArray(xx,&xxa);
1179: BSperm_dvec(xxa,xworka,bsif->pA->perm);CHKERRBS(0);
1180: VecRestoreArray(bsif->xwork,&xworka);
1181: VecRestoreArray(xx,&xxa);
1182: VecPointwiseDivide(xx,bsif->xwork,bsif->diag);
1183: }
1185: VecGetArray(xx,&xxa);
1186: VecGetArray(yy,&yya);
1187: /* Do lower triangular multiplication: [ y = L * xwork ] */
1188: if (bspinfo->single) {
1189: BSforward1(bsif->pA,xxa,yya,bsif->comm_pA,bspinfo);CHKERRBS(0);
1190: } else {
1191: BSforward(bsif->pA,xxa,yya,bsif->comm_pA,bspinfo);CHKERRBS(0);
1192: }
1193:
1194: /* Do upper triangular multiplication: [ y = y + L^{T} * xwork ] */
1195: if (mat->symmetric) {
1196: if (bspinfo->single){
1197: BSbackward1(bsif->pA,xxa,yya,bsif->comm_pA,bspinfo);CHKERRBS(0);
1198: } else {
1199: BSbackward(bsif->pA,xxa,yya,bsif->comm_pA,bspinfo);CHKERRBS(0);
1200: }
1201: }
1202: /* not needed for ILU version since forward does it all */
1203: VecRestoreArray(xx,&xxa);
1204: VecRestoreArray(yy,&yya);
1206: /* Apply diagonal scaling to vector: [ y = D^{1/2} * y ] */
1207: if (!bsif->vecs_permscale) {
1208: VecGetArray(bsif->xwork,&xworka);
1209: VecGetArray(xx,&xxa);
1210: BSiperm_dvec(xworka,xxa,bsif->pA->perm);CHKERRBS(0);
1211: VecRestoreArray(bsif->xwork,&xworka);
1212: VecRestoreArray(xx,&xxa);
1213: VecPointwiseDivide(bsif->xwork,yy,bsif->diag);
1214: VecGetArray(bsif->xwork,&xworka);
1215: VecGetArray(yy,&yya);
1216: BSiperm_dvec(xworka,yya,bsif->pA->perm);CHKERRBS(0);
1217: VecRestoreArray(bsif->xwork,&xworka);
1218: VecRestoreArray(yy,&yya);
1219: }
1220: PetscLogFlops(2*bsif->nz - mat->cmap.n);
1222: return(0);
1223: }
1227: PetscErrorCode MatMultAdd_MPIRowbs(Mat mat,Vec xx,Vec yy,Vec zz)
1228: {
1230: PetscScalar one = 1.0;
1233: (*mat->ops->mult)(mat,xx,zz);
1234: VecAXPY(zz,one,yy);
1235: return(0);
1236: }
1240: PetscErrorCode MatGetInfo_MPIRowbs(Mat A,MatInfoType flag,MatInfo *info)
1241: {
1242: Mat_MPIRowbs *mat = (Mat_MPIRowbs*)A->data;
1243: PetscReal isend[5],irecv[5];
1247: info->rows_global = (double)A->rmap.N;
1248: info->columns_global = (double)A->cmap.N;
1249: info->rows_local = (double)A->cmap.n;
1250: info->columns_local = (double)A->rmap.n;
1251: info->block_size = 1.0;
1252: info->mallocs = (double)mat->reallocs;
1253: isend[0] = mat->nz; isend[1] = mat->maxnz; isend[2] = mat->maxnz - mat->nz;
1254: isend[3] = ((PetscObject)A)->mem; isend[4] = info->mallocs;
1256: if (flag == MAT_LOCAL) {
1257: info->nz_used = isend[0];
1258: info->nz_allocated = isend[1];
1259: info->nz_unneeded = isend[2];
1260: info->memory = isend[3];
1261: info->mallocs = isend[4];
1262: } else if (flag == MAT_GLOBAL_MAX) {
1263: MPI_Allreduce(isend,irecv,3,MPIU_REAL,MPI_MAX,((PetscObject)A)->comm);
1264: info->nz_used = irecv[0];
1265: info->nz_allocated = irecv[1];
1266: info->nz_unneeded = irecv[2];
1267: info->memory = irecv[3];
1268: info->mallocs = irecv[4];
1269: } else if (flag == MAT_GLOBAL_SUM) {
1270: MPI_Allreduce(isend,irecv,3,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);
1271: info->nz_used = irecv[0];
1272: info->nz_allocated = irecv[1];
1273: info->nz_unneeded = irecv[2];
1274: info->memory = irecv[3];
1275: info->mallocs = irecv[4];
1276: }
1277: return(0);
1278: }
1282: PetscErrorCode MatGetDiagonal_MPIRowbs(Mat mat,Vec v)
1283: {
1284: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
1285: BSsprow **rs = a->A->rows;
1287: int i,n;
1288: PetscScalar *x,zero = 0.0;
1291: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1292: if (!a->blocksolveassembly) {
1293: MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
1294: }
1296: VecSet(v,zero);
1297: VecGetLocalSize(v,&n);
1298: if (n != mat->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1299: VecGetArray(v,&x);
1300: for (i=0; i<mat->rmap.n; i++) {
1301: x[i] = rs[i]->nz[rs[i]->diag_ind];
1302: }
1303: VecRestoreArray(v,&x);
1304: return(0);
1305: }
1309: PetscErrorCode MatDestroy_MPIRowbs(Mat mat)
1310: {
1311: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
1312: BSspmat *A = a->A;
1313: BSsprow *vs;
1315: int i;
1318: #if defined(PETSC_USE_LOG)
1319: PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->rmap.N,mat->cmap.N);
1320: #endif
1321: MatStashDestroy_Private(&mat->stash);
1322: if (a->bsmap) {
1323: PetscFree(a->bsmap->vlocal2global);
1324: PetscFree(a->bsmap->vglobal2local);
1325: if (a->bsmap->vglobal2proc) (*a->bsmap->free_g2p)(a->bsmap->vglobal2proc);
1326: PetscFree(a->bsmap);
1327: }
1329: if (A) {
1330: for (i=0; i<mat->rmap.n; i++) {
1331: vs = A->rows[i];
1332: MatFreeRowbs_Private(mat,vs->length,vs->col,vs->nz);
1333: }
1334: /* Note: A->map = a->bsmap is freed above */
1335: PetscFree(A->rows);
1336: PetscFree(A);
1337: }
1338: if (a->procinfo) {BSfree_ctx(a->procinfo);CHKERRBS(0);}
1339: if (a->diag) {VecDestroy(a->diag);}
1340: if (a->xwork) {VecDestroy(a->xwork);}
1341: if (a->pA) {BSfree_par_mat(a->pA);CHKERRBS(0);}
1342: if (a->fpA) {BSfree_copy_par_mat(a->fpA);CHKERRBS(0);}
1343: if (a->comm_pA) {BSfree_comm(a->comm_pA);CHKERRBS(0);}
1344: if (a->comm_fpA) {BSfree_comm(a->comm_fpA);CHKERRBS(0);}
1345: PetscFree(a->imax);
1346: MPI_Comm_free(&(a->comm_mpirowbs));
1347: PetscFree(a);
1349: PetscObjectChangeTypeName((PetscObject)mat,0);
1350: PetscObjectComposeFunction((PetscObject)mat,"MatMPIRowbsSetPreallocation_C","",PETSC_NULL);
1351: return(0);
1352: }
1356: PetscErrorCode MatSetOption_MPIRowbs(Mat A,MatOption op,PetscTruth flg)
1357: {
1358: Mat_MPIRowbs *a = (Mat_MPIRowbs*)A->data;
1362: switch (op) {
1363: case MAT_ROW_ORIENTED:
1364: a->roworiented = flg;
1365: break;
1366: case MAT_NEW_NONZERO_LOCATIONS:
1367: a->nonew = (flg ? 0 : 1);
1368: break;
1369: case MAT_USE_INODES:
1370: a->bs_color_single = (flg ? 0 : 1);
1371: break;
1372: case MAT_NEW_DIAGONALS:
1373: case MAT_NEW_NONZERO_LOCATION_ERR:
1374: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1375: case MAT_USE_HASH_TABLE:
1376: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1377: break;
1378: case MAT_IGNORE_OFF_PROC_ENTRIES:
1379: a->donotstash = flg;
1380: break;
1381: case MAT_KEEP_ZEROED_ROWS:
1382: a->keepzeroedrows = flg;
1383: break;
1384: case MAT_SYMMETRIC:
1385: BSset_mat_symmetric(a->A,PETSC_TRUE);CHKERRBS(0);
1386: break;
1387: case MAT_STRUCTURALLY_SYMMETRIC:
1388: case MAT_HERMITIAN:
1389: case MAT_SYMMETRY_ETERNAL:
1390: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1391: break;
1392: default:
1393: SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
1394: break;
1395: }
1396: return(0);
1397: }
1401: PetscErrorCode MatGetRow_MPIRowbs(Mat AA,int row,int *nz,int **idx,PetscScalar **v)
1402: {
1403: Mat_MPIRowbs *mat = (Mat_MPIRowbs*)AA->data;
1404: BSspmat *A = mat->A;
1405: BSsprow *rs;
1406:
1408: if (row < AA->rmap.rstart || row >= AA->rmap.rend) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Only local rows");
1410: rs = A->rows[row - AA->rmap.rstart];
1411: *nz = rs->length;
1412: if (v) *v = rs->nz;
1413: if (idx) *idx = rs->col;
1414: return(0);
1415: }
1419: PetscErrorCode MatRestoreRow_MPIRowbs(Mat A,int row,int *nz,int **idx,PetscScalar **v)
1420: {
1422: return(0);
1423: }
1425: /* ------------------------------------------------------------------ */
1429: PetscErrorCode MatSetUpPreallocation_MPIRowbs(Mat A)
1430: {
1434: MatMPIRowbsSetPreallocation(A,PETSC_DEFAULT,0);
1435: return(0);
1436: }
1438: /* -------------------------------------------------------------------*/
1439: static struct _MatOps MatOps_Values = {MatSetValues_MPIRowbs,
1440: MatGetRow_MPIRowbs,
1441: MatRestoreRow_MPIRowbs,
1442: MatMult_MPIRowbs,
1443: /* 4*/ MatMultAdd_MPIRowbs,
1444: MatMult_MPIRowbs,
1445: MatMultAdd_MPIRowbs,
1446: MatSolve_MPIRowbs,
1447: 0,
1448: 0,
1449: /*10*/ 0,
1450: 0,
1451: 0,
1452: 0,
1453: 0,
1454: /*15*/ MatGetInfo_MPIRowbs,
1455: 0,
1456: MatGetDiagonal_MPIRowbs,
1457: 0,
1458: MatNorm_MPIRowbs,
1459: /*20*/ MatAssemblyBegin_MPIRowbs,
1460: MatAssemblyEnd_MPIRowbs,
1461: 0,
1462: MatSetOption_MPIRowbs,
1463: MatZeroEntries_MPIRowbs,
1464: /*25*/ MatZeroRows_MPIRowbs,
1465: 0,
1466: MatLUFactorNumeric_MPIRowbs,
1467: 0,
1468: MatCholeskyFactorNumeric_MPIRowbs,
1469: /*30*/ MatSetUpPreallocation_MPIRowbs,
1470: MatILUFactorSymbolic_MPIRowbs,
1471: MatIncompleteCholeskyFactorSymbolic_MPIRowbs,
1472: 0,
1473: 0,
1474: /*35*/ 0,
1475: MatForwardSolve_MPIRowbs,
1476: MatBackwardSolve_MPIRowbs,
1477: 0,
1478: 0,
1479: /*40*/ 0,
1480: MatGetSubMatrices_MPIRowbs,
1481: 0,
1482: 0,
1483: 0,
1484: /*45*/ 0,
1485: MatScale_MPIRowbs,
1486: 0,
1487: 0,
1488: 0,
1489: /*50*/ 0,
1490: 0,
1491: 0,
1492: 0,
1493: 0,
1494: /*55*/ 0,
1495: 0,
1496: 0,
1497: 0,
1498: 0,
1499: /*60*/ MatGetSubMatrix_MPIRowbs,
1500: MatDestroy_MPIRowbs,
1501: MatView_MPIRowbs,
1502: 0,
1503: MatUseScaledForm_MPIRowbs,
1504: /*65*/ MatScaleSystem_MPIRowbs,
1505: MatUnScaleSystem_MPIRowbs,
1506: 0,
1507: 0,
1508: 0,
1509: /*70*/ 0,
1510: 0,
1511: 0,
1512: 0,
1513: 0,
1514: /*75*/ 0,
1515: 0,
1516: 0,
1517: 0,
1518: 0,
1519: /*80*/ 0,
1520: 0,
1521: 0,
1522: 0,
1523: MatLoad_MPIRowbs,
1524: /*85*/ 0,
1525: 0,
1526: 0,
1527: 0,
1528: 0,
1529: /*90*/ 0,
1530: 0,
1531: 0,
1532: 0,
1533: 0,
1534: /*95*/ 0,
1535: 0,
1536: 0,
1537: 0};
1539: /* ------------------------------------------------------------------- */
1544: PetscErrorCode MatMPIRowbsSetPreallocation_MPIRowbs(Mat mat,int nz,const int nnz[])
1545: {
1549: mat->preallocated = PETSC_TRUE;
1550: MatCreateMPIRowbs_local(mat,nz,nnz);
1551: return(0);
1552: }
1555: /*MC
1556: MATMPIROWBS - MATMPIROWBS = "mpirowbs" - A matrix type providing ILU and ICC for distributed sparse matrices for use
1557: with the external package BlockSolve95. If BlockSolve95 is installed (see the manual for instructions
1558: on how to declare the existence of external packages), a matrix type can be constructed which invokes
1559: BlockSolve95 preconditioners and solvers.
1561: Options Database Keys:
1562: . -mat_type mpirowbs - sets the matrix type to "mpirowbs" during a call to MatSetFromOptions()
1564: Level: beginner
1566: .seealso: MatCreateMPIRowbs
1567: M*/
1572: PetscErrorCode MatCreate_MPIRowbs(Mat A)
1573: {
1574: Mat_MPIRowbs *a;
1575: BSmapping *bsmap;
1576: BSoff_map *bsoff;
1578: int *offset,m,M;
1579: PetscTruth flg1,flg3;
1580: BSprocinfo *bspinfo;
1581: MPI_Comm comm;
1582:
1584: comm = ((PetscObject)A)->comm;
1586: PetscNewLog(A,Mat_MPIRowbs,&a);
1587: A->data = (void*)a;
1588: PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));
1589: A->factor = 0;
1590: A->mapping = 0;
1591: a->vecs_permscale = PETSC_FALSE;
1592: A->insertmode = NOT_SET_VALUES;
1593: a->blocksolveassembly = 0;
1594: a->keepzeroedrows = PETSC_FALSE;
1596: MPI_Comm_rank(comm,&a->rank);
1597: MPI_Comm_size(comm,&a->size);
1600: PetscMapSetBlockSize(&A->rmap,1);
1601: PetscMapSetBlockSize(&A->cmap,1);
1602: PetscMapSetUp(&A->rmap);
1603: PetscMapSetUp(&A->cmap);
1604: m = A->rmap.n;
1605: M = A->rmap.N;
1607: PetscMalloc((A->rmap.n+1)*sizeof(int),&a->imax);
1608: a->reallocs = 0;
1610: /* build cache for off array entries formed */
1611: MatStashCreate_Private(((PetscObject)A)->comm,1,&A->stash);
1612: a->donotstash = PETSC_FALSE;
1614: /* Initialize BlockSolve information */
1615: a->A = 0;
1616: a->pA = 0;
1617: a->comm_pA = 0;
1618: a->fpA = 0;
1619: a->comm_fpA = 0;
1620: a->alpha = 1.0;
1621: a->0;
1622: a->failures = 0;
1623: MPI_Comm_dup(((PetscObject)A)->comm,&(a->comm_mpirowbs));
1624: VecCreateMPI(((PetscObject)A)->comm,A->rmap.n,A->rmap.N,&(a->diag));
1625: VecDuplicate(a->diag,&(a->xwork));
1626: PetscLogObjectParent(A,a->diag); PetscLogObjectParent(A,a->xwork);
1627: PetscLogObjectMemory(A,(A->rmap.n+1)*sizeof(PetscScalar));
1628: bspinfo = BScreate_ctx();CHKERRBS(0);
1629: a->procinfo = bspinfo;
1630: BSctx_set_id(bspinfo,a->rank);CHKERRBS(0);
1631: BSctx_set_np(bspinfo,a->size);CHKERRBS(0);
1632: BSctx_set_ps(bspinfo,a->comm_mpirowbs);CHKERRBS(0);
1633: BSctx_set_cs(bspinfo,INT_MAX);CHKERRBS(0);
1634: BSctx_set_is(bspinfo,INT_MAX);CHKERRBS(0);
1635: BSctx_set_ct(bspinfo,IDO);CHKERRBS(0);
1636: #if defined(PETSC_USE_DEBUG)
1637: BSctx_set_err(bspinfo,1);CHKERRBS(0); /* BS error checking */
1638: #endif
1639: BSctx_set_rt(bspinfo,1);CHKERRBS(0);
1640: #if defined (PETSC_USE_INFO)
1641: PetscOptionsHasName(PETSC_NULL,"-info",&flg1);
1642: if (flg1) {
1643: BSctx_set_pr(bspinfo,1);CHKERRBS(0);
1644: }
1645: #endif
1646: PetscOptionsBegin(((PetscObject)A)->comm,PETSC_NULL,"Options for MPIROWBS matrix","Mat");
1647: PetscOptionsTruth("-pc_factor_factorpointwise","Do not optimize for inodes (slow)",PETSC_NULL,PETSC_FALSE,&flg1,PETSC_NULL);
1648: PetscOptionsTruth("-mat_rowbs_no_inode","Do not optimize for inodes (slow)",PETSC_NULL,PETSC_FALSE,&flg3,PETSC_NULL);
1649: PetscOptionsEnd();
1650: if (flg1 || flg3) {
1651: BSctx_set_si(bspinfo,1);CHKERRBS(0);
1652: } else {
1653: BSctx_set_si(bspinfo,0);CHKERRBS(0);
1654: }
1655: #if defined(PETSC_USE_LOG)
1656: MLOG_INIT(); /* Initialize logging */
1657: #endif
1659: /* Compute global offsets */
1660: offset = &A->rmap.rstart;
1662: PetscNewLog(A,BSmapping,&a->bsmap);
1663: bsmap = a->bsmap;
1664: PetscMalloc(sizeof(int),&bsmap->vlocal2global);
1665: *((int*)bsmap->vlocal2global) = (*offset);
1666: bsmap->flocal2global = BSloc2glob;
1667: bsmap->free_l2g = 0;
1668: PetscMalloc(sizeof(int),&bsmap->vglobal2local);
1669: *((int*)bsmap->vglobal2local) = (*offset);
1670: bsmap->fglobal2local = BSglob2loc;
1671: bsmap->free_g2l = 0;
1672: bsoff = BSmake_off_map(*offset,bspinfo,A->rmap.N);
1673: bsmap->vglobal2proc = (void*)bsoff;
1674: bsmap->fglobal2proc = BSglob2proc;
1675: bsmap->free_g2p = (void(*)(void*)) BSfree_off_map;
1676: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMPIRowbsSetPreallocation_C",
1677: "MatMPIRowbsSetPreallocation_MPIRowbs",
1678: MatMPIRowbsSetPreallocation_MPIRowbs);
1679: PetscObjectChangeTypeName((PetscObject)A,MATMPIROWBS);
1680: return(0);
1681: }
1686: /* @
1687: MatMPIRowbsSetPreallocation - Sets the number of expected nonzeros
1688: per row in the matrix.
1690: Input Parameter:
1691: + mat - matrix
1692: . nz - maximum expected for any row
1693: - nzz - number expected in each row
1695: Note:
1696: This routine is valid only for matrices stored in the MATMPIROWBS
1697: format.
1698: @ */
1699: PetscErrorCode MatMPIRowbsSetPreallocation(Mat mat,int nz,const int nnz[])
1700: {
1701: PetscErrorCode ierr,(*f)(Mat,int,const int[]);
1704: PetscObjectQueryFunction((PetscObject)mat,"MatMPIRowbsSetPreallocation_C",(void (**)(void))&f);
1705: if (f) {
1706: (*f)(mat,nz,nnz);
1707: }
1708: return(0);
1709: }
1711: /* --------------- extra BlockSolve-specific routines -------------- */
1714: /* @
1715: MatGetBSProcinfo - Gets the BlockSolve BSprocinfo context, which the
1716: user can then manipulate to alter the default parameters.
1718: Input Parameter:
1719: mat - matrix
1721: Output Parameter:
1722: procinfo - processor information context
1724: Note:
1725: This routine is valid only for matrices stored in the MATMPIROWBS
1726: format.
1727: @ */
1728: PetscErrorCode MatGetBSProcinfo(Mat mat,BSprocinfo *procinfo)
1729: {
1730: Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
1731: PetscTruth ismpirowbs;
1735: PetscTypeCompare((PetscObject)mat,MATMPIROWBS,&ismpirowbs);
1736: if (!ismpirowbs) SETERRQ(PETSC_ERR_ARG_WRONG,"For MATMPIROWBS matrix type");
1737: procinfo = a->procinfo;
1738: return(0);
1739: }
1743: PetscErrorCode MatLoad_MPIRowbs(PetscViewer viewer,MatType type,Mat *newmat)
1744: {
1745: Mat_MPIRowbs *a;
1746: BSspmat *A;
1747: BSsprow **rs;
1748: Mat mat;
1750: int i,nz,j,rstart,rend,fd,*ourlens,*sndcounts = 0,*procsnz;
1751: int header[4],rank,size,*rowlengths = 0,M,m,*rowners,maxnz,*cols;
1752: PetscScalar *vals;
1753: MPI_Comm comm = ((PetscObject)viewer)->comm;
1754: MPI_Status status;
1757: MPI_Comm_size(comm,&size);
1758: MPI_Comm_rank(comm,&rank);
1759: if (!rank) {
1760: PetscViewerBinaryGetDescriptor(viewer,&fd);
1761: PetscBinaryRead(fd,(char *)header,4,PETSC_INT);
1762: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
1763: if (header[3] < 0) {
1764: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format,cannot load as MPIRowbs");
1765: }
1766: }
1768: MPI_Bcast(header+1,3,MPI_INT,0,comm);
1769: M = header[1];
1771: /* determine ownership of all rows */
1772: m = M/size + ((M % size) > rank);
1773: PetscMalloc((size+2)*sizeof(int),&rowners);
1774: MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1775: rowners[0] = 0;
1776: for (i=2; i<=size; i++) {
1777: rowners[i] += rowners[i-1];
1778: }
1779: rstart = rowners[rank];
1780: rend = rowners[rank+1];
1782: /* distribute row lengths to all processors */
1783: PetscMalloc((rend-rstart)*sizeof(int),&ourlens);
1784: if (!rank) {
1785: PetscMalloc(M*sizeof(int),&rowlengths);
1786: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
1787: PetscMalloc(size*sizeof(int),&sndcounts);
1788: for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1789: MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
1790: PetscFree(sndcounts);
1791: } else {
1792: MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
1793: }
1795: /* create our matrix */
1796: MatCreate(comm,newmat);
1797: MatSetSizes(*newmat,m,m,M,M);
1798: MatSetType(*newmat,type);
1799: MatMPIRowbsSetPreallocation(*newmat,0,ourlens);
1800: mat = *newmat;
1801: PetscFree(ourlens);
1803: a = (Mat_MPIRowbs*)mat->data;
1804: A = a->A;
1805: rs = A->rows;
1807: if (!rank) {
1808: /* calculate the number of nonzeros on each processor */
1809: PetscMalloc(size*sizeof(int),&procsnz);
1810: PetscMemzero(procsnz,size*sizeof(int));
1811: for (i=0; i<size; i++) {
1812: for (j=rowners[i]; j< rowners[i+1]; j++) {
1813: procsnz[i] += rowlengths[j];
1814: }
1815: }
1816: PetscFree(rowlengths);
1818: /* determine max buffer needed and allocate it */
1819: maxnz = 0;
1820: for (i=0; i<size; i++) {
1821: maxnz = PetscMax(maxnz,procsnz[i]);
1822: }
1823: PetscMalloc(maxnz*sizeof(int),&cols);
1825: /* read in my part of the matrix column indices */
1826: nz = procsnz[0];
1827: PetscBinaryRead(fd,cols,nz,PETSC_INT);
1828:
1829: /* insert it into my part of matrix */
1830: nz = 0;
1831: for (i=0; i<A->num_rows; i++) {
1832: for (j=0; j<a->imax[i]; j++) {
1833: rs[i]->col[j] = cols[nz++];
1834: }
1835: rs[i]->length = a->imax[i];
1836: }
1837: /* read in parts for all other processors */
1838: for (i=1; i<size; i++) {
1839: nz = procsnz[i];
1840: PetscBinaryRead(fd,cols,nz,PETSC_INT);
1841: MPI_Send(cols,nz,MPI_INT,i,((PetscObject)mat)->tag,comm);
1842: }
1843: PetscFree(cols);
1844: PetscMalloc(maxnz*sizeof(PetscScalar),&vals);
1846: /* read in my part of the matrix numerical values */
1847: nz = procsnz[0];
1848: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1849:
1850: /* insert it into my part of matrix */
1851: nz = 0;
1852: for (i=0; i<A->num_rows; i++) {
1853: for (j=0; j<a->imax[i]; j++) {
1854: rs[i]->nz[j] = vals[nz++];
1855: }
1856: }
1857: /* read in parts for all other processors */
1858: for (i=1; i<size; i++) {
1859: nz = procsnz[i];
1860: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1861: MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)mat)->tag,comm);
1862: }
1863: PetscFree(vals);
1864: PetscFree(procsnz);
1865: } else {
1866: /* determine buffer space needed for message */
1867: nz = 0;
1868: for (i=0; i<A->num_rows; i++) {
1869: nz += a->imax[i];
1870: }
1871: PetscMalloc(nz*sizeof(int),&cols);
1873: /* receive message of column indices*/
1874: MPI_Recv(cols,nz,MPI_INT,0,((PetscObject)mat)->tag,comm,&status);
1875: MPI_Get_count(&status,MPI_INT,&maxnz);
1876: if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong");
1878: /* insert it into my part of matrix */
1879: nz = 0;
1880: for (i=0; i<A->num_rows; i++) {
1881: for (j=0; j<a->imax[i]; j++) {
1882: rs[i]->col[j] = cols[nz++];
1883: }
1884: rs[i]->length = a->imax[i];
1885: }
1886: PetscFree(cols);
1887: PetscMalloc(nz*sizeof(PetscScalar),&vals);
1889: /* receive message of values*/
1890: MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)mat)->tag,comm,&status);
1891: MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
1892: if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong");
1894: /* insert it into my part of matrix */
1895: nz = 0;
1896: for (i=0; i<A->num_rows; i++) {
1897: for (j=0; j<a->imax[i]; j++) {
1898: rs[i]->nz[j] = vals[nz++];
1899: }
1900: rs[i]->length = a->imax[i];
1901: }
1902: PetscFree(vals);
1903: }
1904: PetscFree(rowners);
1905: a->nz = a->maxnz;
1906: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
1907: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
1908: return(0);
1909: }
1911: /*
1912: Special destroy and view routines for factored matrices
1913: */
1916: static PetscErrorCode MatDestroy_MPIRowbs_Factored(Mat mat)
1917: {
1919: #if defined(PETSC_USE_LOG)
1920: PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->rmap.N,mat->cmap.N);
1921: #endif
1922: return(0);
1923: }
1927: static PetscErrorCode MatView_MPIRowbs_Factored(Mat mat,PetscViewer viewer)
1928: {
1932: MatView((Mat) mat->data,viewer);
1933: return(0);
1934: }
1938: PetscErrorCode MatIncompleteCholeskyFactorSymbolic_MPIRowbs(Mat mat,IS isrow,MatFactorInfo *info,Mat *newfact)
1939: {
1940: /* Note: f is not currently used in BlockSolve */
1941: Mat newmat;
1942: Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)mat->data;
1944: PetscTruth idn;
1947: if (isrow) {
1948: ISIdentity(isrow,&idn);
1949: if (!idn) SETERRQ(PETSC_ERR_SUP,"Only identity row permutation supported");
1950: }
1952: if (!mat->symmetric) {
1953: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"To use incomplete Cholesky \n\
1954: preconditioning with a MATMPIROWBS matrix you must declare it to be \n\
1955: symmetric using the option MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)");
1956: }
1958: /* If the icc_storage flag wasn't set before the last blocksolveassembly, */
1959: /* we must completely redo the assembly as a different storage format is required. */
1960: if (mbs->blocksolveassembly && !mbs->assembled_icc_storage) {
1961: mat->same_nonzero = PETSC_FALSE;
1962: mbs->blocksolveassembly = 0;
1963: }
1965: if (!mbs->blocksolveassembly) {
1966: BSset_mat_icc_storage(mbs->A,PETSC_TRUE);CHKERRBS(0);
1967: BSset_mat_symmetric(mbs->A,PETSC_TRUE);CHKERRBS(0);
1968: MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
1969: }
1971: /* Copy permuted matrix */
1972: if (mbs->fpA) {BSfree_copy_par_mat(mbs->fpA);CHKERRBS(0);}
1973: mbs->fpA = BScopy_par_mat(mbs->pA);CHKERRBS(0);
1975: /* Set up the communication for factorization */
1976: if (mbs->comm_fpA) {BSfree_comm(mbs->comm_fpA);CHKERRBS(0);}
1977: mbs->comm_fpA = BSsetup_factor(mbs->fpA,mbs->procinfo);CHKERRBS(0);
1979: /*
1980: Create a new Mat structure to hold the "factored" matrix,
1981: not this merely contains a pointer to the original matrix, since
1982: the original matrix contains the factor information.
1983: */
1984: PetscHeaderCreate(newmat,_p_Mat,struct _MatOps,MAT_COOKIE,-1,"Mat",((PetscObject)mat)->comm,MatDestroy,MatView);
1986: newmat->data = (void*)mat;
1987: PetscMemcpy(newmat->ops,&MatOps_Values,sizeof(struct _MatOps));
1988: newmat->ops->destroy = MatDestroy_MPIRowbs_Factored;
1989: newmat->ops->view = MatView_MPIRowbs_Factored;
1990: newmat->factor = 1;
1991: newmat->preallocated = PETSC_TRUE;
1992: PetscMapCopy(((PetscObject)mat)->comm,&mat->rmap,&newmat->rmap);
1993: PetscMapCopy(((PetscObject)mat)->comm,&mat->cmap,&newmat->cmap);
1995: PetscStrallocpy(MATMPIROWBS,&((PetscObject)newmat)->type_name);
1997: *newfact = newmat;
1998: return(0);
1999: }
2003: PetscErrorCode MatILUFactorSymbolic_MPIRowbs(Mat mat,IS isrow,IS iscol,MatFactorInfo* info,Mat *newfact)
2004: {
2005: Mat newmat;
2006: Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)mat->data;
2008: PetscTruth idn;
2011: if (info->levels) SETERRQ(PETSC_ERR_SUP,"Blocksolve ILU only supports 0 fill");
2012: if (isrow) {
2013: ISIdentity(isrow,&idn);
2014: if (!idn) SETERRQ(PETSC_ERR_SUP,"Only identity row permutation supported");
2015: }
2016: if (iscol) {
2017: ISIdentity(iscol,&idn);
2018: if (!idn) SETERRQ(PETSC_ERR_SUP,"Only identity column permutation supported");
2019: }
2021: if (!mbs->blocksolveassembly) {
2022: MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
2023: }
2024:
2025: /* if (mat->symmetric) { */
2026: /* SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"To use ILU preconditioner with \n\ */
2027: /* MatCreateMPIRowbs() matrix you CANNOT declare it to be a symmetric matrix\n\ */
2028: /* using the option MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)"); */
2029: /* } */
2031: /* Copy permuted matrix */
2032: if (mbs->fpA) {BSfree_copy_par_mat(mbs->fpA);CHKERRBS(0);}
2033: mbs->fpA = BScopy_par_mat(mbs->pA);CHKERRBS(0);
2035: /* Set up the communication for factorization */
2036: if (mbs->comm_fpA) {BSfree_comm(mbs->comm_fpA);CHKERRBS(0);}
2037: mbs->comm_fpA = BSsetup_factor(mbs->fpA,mbs->procinfo);CHKERRBS(0);
2039: /*
2040: Create a new Mat structure to hold the "factored" matrix,
2041: not this merely contains a pointer to the original matrix, since
2042: the original matrix contains the factor information.
2043: */
2044: PetscHeaderCreate(newmat,_p_Mat,struct _MatOps,MAT_COOKIE,-1,"Mat",((PetscObject)mat)->comm,MatDestroy,MatView);
2046: newmat->data = (void*)mat;
2047: PetscMemcpy(newmat->ops,&MatOps_Values,sizeof(struct _MatOps));
2048: newmat->ops->destroy = MatDestroy_MPIRowbs_Factored;
2049: newmat->ops->view = MatView_MPIRowbs_Factored;
2050: newmat->factor = 1;
2051: newmat->preallocated = PETSC_TRUE;
2053: PetscMapCopy(((PetscObject)mat)->comm,&mat->rmap,&newmat->rmap);
2054: PetscMapCopy(((PetscObject)mat)->comm,&mat->cmap,&newmat->cmap);
2056: PetscStrallocpy(MATMPIROWBS,&((PetscObject)newmat)->type_name);
2058: *newfact = newmat;
2059: return(0);
2060: }
2064: /*@C
2065: MatCreateMPIRowbs - Creates a sparse parallel matrix in the MATMPIROWBS
2066: format. This format is intended primarily as an interface for BlockSolve95.
2068: Collective on MPI_Comm
2070: Input Parameters:
2071: + comm - MPI communicator
2072: . m - number of local rows (or PETSC_DECIDE to have calculated)
2073: . M - number of global rows (or PETSC_DECIDE to have calculated)
2074: . nz - number of nonzeros per row (same for all local rows)
2075: - nnz - number of nonzeros per row (possibly different for each row).
2077: Output Parameter:
2078: . newA - the matrix
2080: Notes:
2081: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2082: than it must be used on all processors that share the object for that argument.
2084: The user MUST specify either the local or global matrix dimensions
2085: (possibly both).
2087: Specify the preallocated storage with either nz or nnz (not both). Set
2088: nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2089: allocation.
2091: Notes:
2092: By default, the matrix is assumed to be nonsymmetric; the user can
2093: take advantage of special optimizations for symmetric matrices by calling
2094: $ MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)
2095: $ MatSetOption(mat,MAT_SYMMETRY_ETERNAL,PETSC_TRUE)
2096: BEFORE calling the routine MatAssemblyBegin().
2098: Internally, the MATMPIROWBS format inserts zero elements to the
2099: matrix if necessary, so that nonsymmetric matrices are considered
2100: to be symmetric in terms of their sparsity structure; this format
2101: is required for use of the parallel communication routines within
2102: BlockSolve95. In particular, if the matrix element A[i,j] exists,
2103: then PETSc will internally allocate a 0 value for the element
2104: A[j,i] during MatAssemblyEnd() if the user has not already set
2105: a value for the matrix element A[j,i].
2107: Options Database Keys:
2108: . -mat_rowbs_no_inode - Do not use inodes.
2110: Level: intermediate
2111:
2112: .keywords: matrix, row, symmetric, sparse, parallel, BlockSolve
2114: .seealso: MatCreate(), MatSetValues()
2115: @*/
2116: PetscErrorCode MatCreateMPIRowbs(MPI_Comm comm,int m,int M,int nz,const int nnz[],Mat *newA)
2117: {
2119:
2121: MatCreate(comm,newA);
2122: MatSetSizes(*newA,m,m,M,M);
2123: MatSetType(*newA,MATMPIROWBS);
2124: MatMPIRowbsSetPreallocation(*newA,nz,nnz);
2125: return(0);
2126: }
2129: /* -------------------------------------------------------------------------*/
2131: #include src/mat/impls/aij/seq/aij.h
2132: #include src/mat/impls/aij/mpi/mpiaij.h
2136: PetscErrorCode MatGetSubMatrices_MPIRowbs(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
2137: {
2139: int nmax,nstages_local,nstages,i,pos,max_no;
2143: /* Allocate memory to hold all the submatrices */
2144: if (scall != MAT_REUSE_MATRIX) {
2145: PetscMalloc((ismax+1)*sizeof(Mat),submat);
2146: }
2147:
2148: /* Determine the number of stages through which submatrices are done */
2149: nmax = 20*1000000 / (C->cmap.N * sizeof(int));
2150: if (!nmax) nmax = 1;
2151: nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
2153: /* Make sure every processor loops through the nstages */
2154: MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,((PetscObject)C)->comm);
2156: for (i=0,pos=0; i<nstages; i++) {
2157: if (pos+nmax <= ismax) max_no = nmax;
2158: else if (pos == ismax) max_no = 0;
2159: else max_no = ismax-pos;
2160: MatGetSubMatrices_MPIRowbs_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);
2161: pos += max_no;
2162: }
2163: return(0);
2164: }
2165: /* -------------------------------------------------------------------------*/
2166: /* for now MatGetSubMatrices_MPIRowbs_Local get MPIAij submatrices of input
2167: matrix and preservs zeroes from structural symetry
2168: */
2171: PetscErrorCode MatGetSubMatrices_MPIRowbs_Local(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
2172: {
2173: Mat_MPIRowbs *c = (Mat_MPIRowbs *)(C->data);
2174: BSspmat *A = c->A;
2175: Mat_SeqAIJ *mat;
2177: int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size;
2178: int **sbuf1,**sbuf2,rank,m,i,j,k,l,ct1,ct2,**rbuf1,row,proc;
2179: int nrqs,msz,**ptr,idx,*req_size,*ctr,*pa,*tmp,tcol,nrqr;
2180: int **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2,**rmap;
2181: int **cmap,**lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i;
2182: int len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*cmap_i,*lens_i;
2183: int *rmap_i,tag0,tag1,tag2,tag3;
2184: MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
2185: MPI_Request *r_waits4,*s_waits3,*s_waits4;
2186: MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
2187: MPI_Status *r_status3,*r_status4,*s_status4;
2188: MPI_Comm comm;
2189: FLOAT **rbuf4,**sbuf_aa,*vals,*sbuf_aa_i;
2190: PetscScalar *mat_a;
2191: PetscTruth sorted;
2192: int *onodes1,*olengths1;
2195: comm = ((PetscObject)C)->comm;
2196: tag0 = ((PetscObject)C)->tag;
2197: size = c->size;
2198: rank = c->rank;
2199: m = C->rmap.N;
2200:
2201: /* Get some new tags to keep the communication clean */
2202: PetscObjectGetNewTag((PetscObject)C,&tag1);
2203: PetscObjectGetNewTag((PetscObject)C,&tag2);
2204: PetscObjectGetNewTag((PetscObject)C,&tag3);
2206: /* Check if the col indices are sorted */
2207: for (i=0; i<ismax; i++) {
2208: ISSorted(isrow[i],&sorted);
2209: if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
2210: ISSorted(iscol[i],&sorted);
2211: /* if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted"); */
2212: }
2214: len = (2*ismax+1)*(sizeof(int*)+ sizeof(int)) + (m+1)*sizeof(int);
2215: PetscMalloc(len,&irow);
2216: icol = irow + ismax;
2217: nrow = (int*)(icol + ismax);
2218: ncol = nrow + ismax;
2219: rtable = ncol + ismax;
2221: for (i=0; i<ismax; i++) {
2222: ISGetIndices(isrow[i],&irow[i]);
2223: ISGetIndices(iscol[i],&icol[i]);
2224: ISGetLocalSize(isrow[i],&nrow[i]);
2225: ISGetLocalSize(iscol[i],&ncol[i]);
2226: }
2228: /* Create hash table for the mapping :row -> proc*/
2229: for (i=0,j=0; i<size; i++) {
2230: jmax = C->rmap.range[i+1];
2231: for (; j<jmax; j++) {
2232: rtable[j] = i;
2233: }
2234: }
2236: /* evaluate communication - mesg to who, length of mesg, and buffer space
2237: required. Based on this, buffers are allocated, and data copied into them*/
2238: PetscMalloc(size*4*sizeof(int),&w1); /* mesg size */
2239: w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/
2240: w3 = w2 + size; /* no of IS that needs to be sent to proc i */
2241: w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */
2242: PetscMemzero(w1,size*3*sizeof(int)); /* initialize work vector*/
2243: for (i=0; i<ismax; i++) {
2244: PetscMemzero(w4,size*sizeof(int)); /* initialize work vector*/
2245: jmax = nrow[i];
2246: irow_i = irow[i];
2247: for (j=0; j<jmax; j++) {
2248: row = irow_i[j];
2249: proc = rtable[row];
2250: w4[proc]++;
2251: }
2252: for (j=0; j<size; j++) {
2253: if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
2254: }
2255: }
2256:
2257: nrqs = 0; /* no of outgoing messages */
2258: msz = 0; /* total mesg length (for all procs) */
2259: w1[rank] = 0; /* no mesg sent to self */
2260: w3[rank] = 0;
2261: for (i=0; i<size; i++) {
2262: if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
2263: }
2264: PetscMalloc((nrqs+1)*sizeof(int),&pa); /*(proc -array)*/
2265: for (i=0,j=0; i<size; i++) {
2266: if (w1[i]) { pa[j] = i; j++; }
2267: }
2269: /* Each message would have a header = 1 + 2*(no of IS) + data */
2270: for (i=0; i<nrqs; i++) {
2271: j = pa[i];
2272: w1[j] += w2[j] + 2* w3[j];
2273: msz += w1[j];
2274: }
2276: /* Determine the number of messages to expect, their lengths, from from-ids */
2277: PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);
2278: PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);
2280: /* Now post the Irecvs corresponding to these messages */
2281: PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);
2282:
2283: PetscFree(onodes1);
2284: PetscFree(olengths1);
2285:
2286: /* Allocate Memory for outgoing messages */
2287: len = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int);
2288: PetscMalloc(len,&sbuf1);
2289: ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */
2290: PetscMemzero(sbuf1,2*size*sizeof(int*));
2291: /* allocate memory for outgoing data + buf to receive the first reply */
2292: tmp = (int*)(ptr + size);
2293: ctr = tmp + 2*msz;
2295: {
2296: int *iptr = tmp,ict = 0;
2297: for (i=0; i<nrqs; i++) {
2298: j = pa[i];
2299: iptr += ict;
2300: sbuf1[j] = iptr;
2301: ict = w1[j];
2302: }
2303: }
2305: /* Form the outgoing messages */
2306: /* Initialize the header space */
2307: for (i=0; i<nrqs; i++) {
2308: j = pa[i];
2309: sbuf1[j][0] = 0;
2310: PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));
2311: ptr[j] = sbuf1[j] + 2*w3[j] + 1;
2312: }
2313:
2314: /* Parse the isrow and copy data into outbuf */
2315: for (i=0; i<ismax; i++) {
2316: PetscMemzero(ctr,size*sizeof(int));
2317: irow_i = irow[i];
2318: jmax = nrow[i];
2319: for (j=0; j<jmax; j++) { /* parse the indices of each IS */
2320: row = irow_i[j];
2321: proc = rtable[row];
2322: if (proc != rank) { /* copy to the outgoing buf*/
2323: ctr[proc]++;
2324: *ptr[proc] = row;
2325: ptr[proc]++;
2326: }
2327: }
2328: /* Update the headers for the current IS */
2329: for (j=0; j<size; j++) { /* Can Optimise this loop too */
2330: if ((ctr_j = ctr[j])) {
2331: sbuf1_j = sbuf1[j];
2332: k = ++sbuf1_j[0];
2333: sbuf1_j[2*k] = ctr_j;
2334: sbuf1_j[2*k-1] = i;
2335: }
2336: }
2337: }
2339: /* Now post the sends */
2340: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
2341: for (i=0; i<nrqs; ++i) {
2342: j = pa[i];
2343: MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);
2344: }
2346: /* Post Receives to capture the buffer size */
2347: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);
2348: PetscMalloc((nrqs+1)*sizeof(int*),&rbuf2);
2349: rbuf2[0] = tmp + msz;
2350: for (i=1; i<nrqs; ++i) {
2351: rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
2352: }
2353: for (i=0; i<nrqs; ++i) {
2354: j = pa[i];
2355: MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);
2356: }
2358: /* Send to other procs the buf size they should allocate */
2359:
2361: /* Receive messages*/
2362: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);
2363: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);
2364: len = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*);
2365: PetscMalloc(len,&sbuf2);
2366: req_size = (int*)(sbuf2 + nrqr);
2367: req_source = req_size + nrqr;
2368:
2369: {
2370: BSsprow **sAi = A->rows;
2371: int id,rstart = C->rmap.rstart;
2372: int *sbuf2_i;
2374: for (i=0; i<nrqr; ++i) {
2375: MPI_Waitany(nrqr,r_waits1,&idx,r_status1+i);
2376: req_size[idx] = 0;
2377: rbuf1_i = rbuf1[idx];
2378: start = 2*rbuf1_i[0] + 1;
2379: MPI_Get_count(r_status1+i,MPI_INT,&end);
2380: PetscMalloc((end+1)*sizeof(int),&sbuf2[idx]);
2381: sbuf2_i = sbuf2[idx];
2382: for (j=start; j<end; j++) {
2383: id = rbuf1_i[j] - rstart;
2384: ncols = (sAi[id])->length;
2385: sbuf2_i[j] = ncols;
2386: req_size[idx] += ncols;
2387: }
2388: req_source[idx] = r_status1[i].MPI_SOURCE;
2389: /* form the header */
2390: sbuf2_i[0] = req_size[idx];
2391: for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
2392: MPI_Isend(sbuf2_i,end,MPI_INT,req_source[idx],tag1,comm,s_waits2+i);
2393: }
2394: }
2395: PetscFree(r_status1);
2396: PetscFree(r_waits1);
2398: /* recv buffer sizes */
2399: /* Receive messages*/
2400:
2401: PetscMalloc((nrqs+1)*sizeof(int*),&rbuf3);
2402: PetscMalloc((nrqs+1)*sizeof(FLOAT *),&rbuf4);
2403: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);
2404: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);
2405: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);
2407: for (i=0; i<nrqs; ++i) {
2408: MPI_Waitany(nrqs,r_waits2,&idx,r_status2+i);
2409: PetscMalloc((rbuf2[idx][0]+1)*sizeof(int),&rbuf3[idx]);
2410: PetscMalloc((rbuf2[idx][0]+1)*sizeof(FLOAT),&rbuf4[idx]);
2411: MPI_Irecv(rbuf3[idx],rbuf2[idx][0],MPI_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idx);
2412: MPI_Irecv(rbuf4[idx],rbuf2[idx][0],MPIU_SCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idx);
2413: }
2414: PetscFree(r_status2);
2415: PetscFree(r_waits2);
2416:
2417: /* Wait on sends1 and sends2 */
2418: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);
2419: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);
2421: if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status1);}
2422: if (nrqr) {MPI_Waitall(nrqr,s_waits2,s_status2);}
2423: PetscFree(s_status1);
2424: PetscFree(s_status2);
2425: PetscFree(s_waits1);
2426: PetscFree(s_waits2);
2428: /* Now allocate buffers for a->j, and send them off */
2429: PetscMalloc((nrqr+1)*sizeof(int*),&sbuf_aj);
2430: for (i=0,j=0; i<nrqr; i++) j += req_size[i];
2431: PetscMalloc((j+1)*sizeof(int),&sbuf_aj[0]);
2432: for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
2433:
2434: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);
2435: {
2436: BSsprow *brow;
2437: int *Acol;
2438: int rstart = C->rmap.rstart;
2440: for (i=0; i<nrqr; i++) {
2441: rbuf1_i = rbuf1[i];
2442: sbuf_aj_i = sbuf_aj[i];
2443: ct1 = 2*rbuf1_i[0] + 1;
2444: ct2 = 0;
2445: for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
2446: kmax = rbuf1[i][2*j];
2447: for (k=0; k<kmax; k++,ct1++) {
2448: brow = A->rows[rbuf1_i[ct1] - rstart];
2449: ncols = brow->length;
2450: Acol = brow->col;
2451: /* load the column indices for this row into cols*/
2452: cols = sbuf_aj_i + ct2;
2453: PetscMemcpy(cols,Acol,ncols*sizeof(int));
2454: /*for (l=0; l<ncols;l++) cols[l]=Acol[l]; */ /* How is it with
2455: mappings?? */
2456: ct2 += ncols;
2457: }
2458: }
2459: MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);
2460: }
2461: }
2462: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);
2463: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);
2465: /* Allocate buffers for a->a, and send them off */
2466: PetscMalloc((nrqr+1)*sizeof(FLOAT*),&sbuf_aa);
2467: for (i=0,j=0; i<nrqr; i++) j += req_size[i];
2468: PetscMalloc((j+1)*sizeof(FLOAT),&sbuf_aa[0]);
2469: for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1];
2470:
2471: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);
2472: {
2473: BSsprow *brow;
2474: FLOAT *Aval;
2475: int rstart = C->rmap.rstart;
2476:
2477: for (i=0; i<nrqr; i++) {
2478: rbuf1_i = rbuf1[i];
2479: sbuf_aa_i = sbuf_aa[i];
2480: ct1 = 2*rbuf1_i[0]+1;
2481: ct2 = 0;
2482: for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
2483: kmax = rbuf1_i[2*j];
2484: for (k=0; k<kmax; k++,ct1++) {
2485: brow = A->rows[rbuf1_i[ct1] - rstart];
2486: ncols = brow->length;
2487: Aval = brow->nz;
2488: /* load the column values for this row into vals*/
2489: vals = sbuf_aa_i+ct2;
2490: PetscMemcpy(vals,Aval,ncols*sizeof(FLOAT));
2491: ct2 += ncols;
2492: }
2493: }
2494: MPI_Isend(sbuf_aa_i,req_size[i],MPIU_SCALAR,req_source[i],tag3,comm,s_waits4+i);
2495: }
2496: }
2497: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);
2498: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);
2499: PetscFree(rbuf1);
2501: /* Form the matrix */
2502: /* create col map */
2503: {
2504: int *icol_i;
2505:
2506: len = (1+ismax)*sizeof(int*)+ ismax*C->cmap.N*sizeof(int);
2507: PetscMalloc(len,&cmap);
2508: cmap[0] = (int*)(cmap + ismax);
2509: PetscMemzero(cmap[0],(1+ismax*C->cmap.N)*sizeof(int));
2510: for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + C->cmap.N; }
2511: for (i=0; i<ismax; i++) {
2512: jmax = ncol[i];
2513: icol_i = icol[i];
2514: cmap_i = cmap[i];
2515: for (j=0; j<jmax; j++) {
2516: cmap_i[icol_i[j]] = j+1;
2517: }
2518: }
2519: }
2521: /* Create lens which is required for MatCreate... */
2522: for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
2523: len = (1+ismax)*sizeof(int*)+ j*sizeof(int);
2524: PetscMalloc(len,&lens);
2525: lens[0] = (int*)(lens + ismax);
2526: PetscMemzero(lens[0],j*sizeof(int));
2527: for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
2528:
2529: /* Update lens from local data */
2530: { BSsprow *Arow;
2531: for (i=0; i<ismax; i++) {
2532: jmax = nrow[i];
2533: cmap_i = cmap[i];
2534: irow_i = irow[i];
2535: lens_i = lens[i];
2536: for (j=0; j<jmax; j++) {
2537: row = irow_i[j];
2538: proc = rtable[row];
2539: if (proc == rank) {
2540: Arow=A->rows[row-C->rmap.rstart];
2541: ncols=Arow->length;
2542: cols=Arow->col;
2543: for (k=0; k<ncols; k++) {
2544: if (cmap_i[cols[k]]) { lens_i[j]++;}
2545: }
2546: }
2547: }
2548: }
2549: }
2550:
2551: /* Create row map*/
2552: len = (1+ismax)*sizeof(int*)+ ismax*C->rmap.N*sizeof(int);
2553: PetscMalloc(len,&rmap);
2554: rmap[0] = (int*)(rmap + ismax);
2555: PetscMemzero(rmap[0],ismax*C->rmap.N*sizeof(int));
2556: for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->rmap.N;}
2557: for (i=0; i<ismax; i++) {
2558: rmap_i = rmap[i];
2559: irow_i = irow[i];
2560: jmax = nrow[i];
2561: for (j=0; j<jmax; j++) {
2562: rmap_i[irow_i[j]] = j;
2563: }
2564: }
2565:
2566: /* Update lens from offproc data */
2567: {
2568: int *rbuf2_i,*rbuf3_i,*sbuf1_i;
2570: for (tmp2=0; tmp2<nrqs; tmp2++) {
2571: MPI_Waitany(nrqs,r_waits3,&i,r_status3+tmp2);
2572: idx = pa[i];
2573: sbuf1_i = sbuf1[idx];
2574: jmax = sbuf1_i[0];
2575: ct1 = 2*jmax+1;
2576: ct2 = 0;
2577: rbuf2_i = rbuf2[i];
2578: rbuf3_i = rbuf3[i];
2579: for (j=1; j<=jmax; j++) {
2580: is_no = sbuf1_i[2*j-1];
2581: max1 = sbuf1_i[2*j];
2582: lens_i = lens[is_no];
2583: cmap_i = cmap[is_no];
2584: rmap_i = rmap[is_no];
2585: for (k=0; k<max1; k++,ct1++) {
2586: row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
2587: max2 = rbuf2_i[ct1];
2588: for (l=0; l<max2; l++,ct2++) {
2589: if (cmap_i[rbuf3_i[ct2]]) {
2590: lens_i[row]++;
2591: }
2592: }
2593: }
2594: }
2595: }
2596: }
2597: PetscFree(r_status3);
2598: PetscFree(r_waits3);
2599: if (nrqr) {MPI_Waitall(nrqr,s_waits3,s_status3);}
2600: PetscFree(s_status3);
2601: PetscFree(s_waits3);
2603: /* Create the submatrices */
2604: if (scall == MAT_REUSE_MATRIX) {
2605: PetscTruth same;
2606:
2607: /*
2608: Assumes new rows are same length as the old rows,hence bug!
2609: */
2610: for (i=0; i<ismax; i++) {
2611: PetscTypeCompare((PetscObject)(submats[i]),MATSEQAIJ,&same);
2612: if (!same) {
2613: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong type");
2614: }
2615: mat = (Mat_SeqAIJ*)(submats[i]->data);
2616: if ((submats[i]->rmap.n != nrow[i]) || (submats[i]->cmap.n != ncol[i])) {
2617: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
2618: }
2619: PetscMemcmp(mat->ilen,lens[i],submats[i]->rmap.n*sizeof(int),&same);
2620: if (!same) {
2621: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
2622: }
2623: /* Initial matrix as if empty */
2624: PetscMemzero(mat->ilen,submats[i]->rmap.n*sizeof(int));
2625: submats[i]->factor = C->factor;
2626: }
2627: } else {
2628: for (i=0; i<ismax; i++) {
2629: /* Here we want to explicitly generate SeqAIJ matrices */
2630: MatCreate(PETSC_COMM_SELF,submats+i);
2631: MatSetSizes(submats[i],nrow[i],ncol[i],nrow[i],ncol[i]);
2632: MatSetType(submats[i],MATSEQAIJ);
2633: MatSeqAIJSetPreallocation(submats[i],0,lens[i]);
2634: }
2635: }
2637: /* Assemble the matrices */
2638: /* First assemble the local rows */
2639: {
2640: int ilen_row,*imat_ilen,*imat_j,*imat_i,old_row;
2641: PetscScalar *imat_a;
2642: BSsprow *Arow;
2643:
2644: for (i=0; i<ismax; i++) {
2645: mat = (Mat_SeqAIJ*)submats[i]->data;
2646: imat_ilen = mat->ilen;
2647: imat_j = mat->j;
2648: imat_i = mat->i;
2649: imat_a = mat->a;
2650: cmap_i = cmap[i];
2651: rmap_i = rmap[i];
2652: irow_i = irow[i];
2653: jmax = nrow[i];
2654: for (j=0; j<jmax; j++) {
2655: row = irow_i[j];
2656: proc = rtable[row];
2657: if (proc == rank) {
2658: old_row = row;
2659: row = rmap_i[row];
2660: ilen_row = imat_ilen[row];
2661:
2662: Arow=A->rows[old_row-C->rmap.rstart];
2663: ncols=Arow->length;
2664: cols=Arow->col;
2665: vals=Arow->nz;
2666:
2667: mat_i = imat_i[row];
2668: mat_a = imat_a + mat_i;
2669: mat_j = imat_j + mat_i;
2670: for (k=0; k<ncols; k++) {
2671: if ((tcol = cmap_i[cols[k]])) {
2672: *mat_j++ = tcol - 1;
2673: *mat_a++ = (PetscScalar)vals[k];
2674: ilen_row++;
2675: }
2676: }
2677: imat_ilen[row] = ilen_row;
2678: }
2679: }
2680: }
2681: }
2683: /* Now assemble the off proc rows*/
2684: {
2685: int *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
2686: int *imat_j,*imat_i;
2687: PetscScalar *imat_a;
2688: FLOAT *rbuf4_i;
2689:
2690: for (tmp2=0; tmp2<nrqs; tmp2++) {
2691: MPI_Waitany(nrqs,r_waits4,&i,r_status4+tmp2);
2692: idx = pa[i];
2693: sbuf1_i = sbuf1[idx];
2694: jmax = sbuf1_i[0];
2695: ct1 = 2*jmax + 1;
2696: ct2 = 0;
2697: rbuf2_i = rbuf2[i];
2698: rbuf3_i = rbuf3[i];
2699: rbuf4_i = rbuf4[i];
2700: for (j=1; j<=jmax; j++) {
2701: is_no = sbuf1_i[2*j-1];
2702: rmap_i = rmap[is_no];
2703: cmap_i = cmap[is_no];
2704: mat = (Mat_SeqAIJ*)submats[is_no]->data;
2705: imat_ilen = mat->ilen;
2706: imat_j = mat->j;
2707: imat_i = mat->i;
2708: imat_a = mat->a;
2709: max1 = sbuf1_i[2*j];
2710: for (k=0; k<max1; k++,ct1++) {
2711: row = sbuf1_i[ct1];
2712: row = rmap_i[row];
2713: ilen = imat_ilen[row];
2714: mat_i = imat_i[row];
2715: mat_a = imat_a + mat_i;
2716: mat_j = imat_j + mat_i;
2717: max2 = rbuf2_i[ct1];
2718: for (l=0; l<max2; l++,ct2++) {
2719: if ((tcol = cmap_i[rbuf3_i[ct2]])) {
2720: *mat_j++ = tcol - 1;
2721: *mat_a++ = (PetscScalar)rbuf4_i[ct2];
2722: ilen++;
2723: }
2724: }
2725: imat_ilen[row] = ilen;
2726: }
2727: }
2728: }
2729: }
2730: PetscFree(r_status4);
2731: PetscFree(r_waits4);
2732: if (nrqr) {MPI_Waitall(nrqr,s_waits4,s_status4);}
2733: PetscFree(s_waits4);
2734: PetscFree(s_status4);
2736: /* Restore the indices */
2737: for (i=0; i<ismax; i++) {
2738: ISRestoreIndices(isrow[i],irow+i);
2739: ISRestoreIndices(iscol[i],icol+i);
2740: }
2742: /* Destroy allocated memory */
2743: PetscFree(irow);
2744: PetscFree(w1);
2745: PetscFree(pa);
2747: PetscFree(sbuf1);
2748: PetscFree(rbuf2);
2749: for (i=0; i<nrqr; ++i) {
2750: PetscFree(sbuf2[i]);
2751: }
2752: for (i=0; i<nrqs; ++i) {
2753: PetscFree(rbuf3[i]);
2754: PetscFree(rbuf4[i]);
2755: }
2757: PetscFree(sbuf2);
2758: PetscFree(rbuf3);
2759: PetscFree(rbuf4);
2760: PetscFree(sbuf_aj[0]);
2761: PetscFree(sbuf_aj);
2762: PetscFree(sbuf_aa[0]);
2763: PetscFree(sbuf_aa);
2764:
2765: PetscFree(cmap);
2766: PetscFree(rmap);
2767: PetscFree(lens);
2769: for (i=0; i<ismax; i++) {
2770: MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);
2771: MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);
2772: }
2773: return(0);
2774: }
2776: /*
2777: can be optimized by send only non-zeroes in iscol IS -
2778: so prebuild submatrix on sending side including A,B partitioning
2779: */
2782: #include src/vec/is/impls/general/general.h
2783: PetscErrorCode MatGetSubMatrix_MPIRowbs(Mat C,IS isrow,IS iscol,int csize,MatReuse scall,Mat *submat)
2784: {
2785: Mat_MPIRowbs *c = (Mat_MPIRowbs*)C->data;
2786: BSspmat *A = c->A;
2787: BSsprow *Arow;
2788: Mat_SeqAIJ *matA,*matB; /* on prac , off proc part of submat */
2789: Mat_MPIAIJ *mat; /* submat->data */
2791: int *irow,*icol,nrow,ncol,*rtable,size,rank,tag0,tag1,tag2,tag3;
2792: int *w1,*w2,*pa,nrqs,nrqr,msz,row_t;
2793: int i,j,k,l,len,jmax,proc,idx;
2794: int **sbuf1,**sbuf2,**rbuf1,**rbuf2,*req_size,**sbuf3,**rbuf3;
2795: FLOAT **rbuf4,**sbuf4; /* FLOAT is from Block Solve 95 library */
2797: int *cmap,*rmap,nlocal,*o_nz,*d_nz,cstart,cend;
2798: int *req_source;
2799: int ncols_t;
2800:
2801:
2802: MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
2803: MPI_Request *r_waits4,*s_waits3,*s_waits4;
2804:
2805: MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
2806: MPI_Status *r_status3,*r_status4,*s_status4;
2807: MPI_Comm comm;
2811: comm = ((PetscObject)C)->comm;
2812: tag0 = ((PetscObject)C)->tag;
2813: size = c->size;
2814: rank = c->rank;
2816: if (size==1) {
2817: if (scall == MAT_REUSE_MATRIX) {
2818: ierr=MatGetSubMatrices(C,1,&isrow,&iscol,MAT_REUSE_MATRIX,&submat);
2819: return(0);
2820: } else {
2821: Mat *newsubmat;
2822:
2823: ierr=MatGetSubMatrices(C,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&newsubmat);
2824: *submat=*newsubmat;
2825: ierr=PetscFree(newsubmat);
2826: return(0);
2827: }
2828: }
2829:
2830: /* Get some new tags to keep the communication clean */
2831: PetscObjectGetNewTag((PetscObject)C,&tag1);
2832: PetscObjectGetNewTag((PetscObject)C,&tag2);
2833: PetscObjectGetNewTag((PetscObject)C,&tag3);
2835: /* Check if the col indices are sorted */
2836: {PetscTruth sorted;
2837: ISSorted(isrow,&sorted);
2838: if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
2839: ISSorted(iscol,&sorted);
2840: if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
2841: }
2842:
2843: ISGetIndices(isrow,&irow);
2844: ISGetIndices(iscol,&icol);
2845: ISGetLocalSize(isrow,&nrow);
2846: ISGetLocalSize(iscol,&ncol);
2847:
2848: if (!isrow) SETERRQ(PETSC_ERR_ARG_SIZ,"Empty ISrow");
2849: if (!iscol) SETERRQ(PETSC_ERR_ARG_SIZ,"Empty IScol");
2850:
2851:
2852: len = (C->rmap.N+1)*sizeof(int);
2853: PetscMalloc(len,&rtable);
2854: /* Create hash table for the mapping :row -> proc*/
2855: for (i=0,j=0; i<size; i++) {
2856: jmax = C->rmap.range[i+1];
2857: for (; j<jmax; j++) {
2858: rtable[j] = i;
2859: }
2860: }
2862: /* evaluate communication - mesg to who, length of mesg, and buffer space
2863: required. Based on this, buffers are allocated, and data copied into them*/
2864: PetscMalloc(size*2*sizeof(int),&w1); /* mesg size */
2865: w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/
2866: PetscMemzero(w1,size*2*sizeof(int)); /* initialize work vector*/
2867: for (j=0; j<nrow; j++) {
2868: row_t = irow[j];
2869: proc = rtable[row_t];
2870: w1[proc]++;
2871: }
2872: nrqs = 0; /* no of outgoing messages */
2873: msz = 0; /* total mesg length (for all procs) */
2874: w1[rank] = 0; /* no mesg sent to self */
2875: for (i=0; i<size; i++) {
2876: if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
2877: }
2878:
2879: PetscMalloc((nrqs+1)*sizeof(int),&pa); /*(proc -array)*/
2880: for (i=0,j=0; i<size; i++) {
2881: if (w1[i]) {
2882: pa[j++] = i;
2883: w1[i]++; /* header for return data */
2884: msz+=w1[i];
2885: }
2886: }
2887:
2888: {int *onodes1,*olengths1;
2889: /* Determine the number of messages to expect, their lengths, from from-ids */
2890: PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);
2891: PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);
2892: /* Now post the Irecvs corresponding to these messages */
2893: PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);
2894: PetscFree(onodes1);
2895: PetscFree(olengths1);
2896: }
2897:
2898: { int **ptr,*iptr,*tmp;
2899: /* Allocate Memory for outgoing messages */
2900: len = 2*size*sizeof(int*) + msz*sizeof(int);
2901: PetscMalloc(len,&sbuf1);
2902: ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */
2903: PetscMemzero(sbuf1,2*size*sizeof(int*));
2904: /* allocate memory for outgoing data + buf to receive the first reply */
2905: tmp = (int*)(ptr + size);
2907: for (i=0,iptr=tmp; i<nrqs; i++) {
2908: j = pa[i];
2909: sbuf1[j] = iptr;
2910: iptr += w1[j];
2911: }
2913: /* Form the outgoing messages */
2914: for (i=0; i<nrqs; i++) {
2915: j = pa[i];
2916: sbuf1[j][0] = 0; /*header */
2917: ptr[j] = sbuf1[j] + 1;
2918: }
2919:
2920: /* Parse the isrow and copy data into outbuf */
2921: for (j=0; j<nrow; j++) {
2922: row_t = irow[j];
2923: proc = rtable[row_t];
2924: if (proc != rank) { /* copy to the outgoing buf*/
2925: sbuf1[proc][0]++;
2926: *ptr[proc] = row_t;
2927: ptr[proc]++;
2928: }
2929: }
2930: } /* block */
2932: /* Now post the sends */
2933:
2934: /* structure of sbuf1[i]/rbuf1[i] : 1 (num of rows) + nrow-local rows (nuberes
2935: * of requested rows)*/
2937: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
2938: for (i=0; i<nrqs; ++i) {
2939: j = pa[i];
2940: MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);
2941: }
2943: /* Post Receives to capture the buffer size */
2944: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);
2945: PetscMalloc((nrqs+1)*sizeof(int*),&rbuf2);
2946: PetscMalloc(msz*sizeof(int)+1,&(rbuf2[0]));
2947: for (i=1; i<nrqs; ++i) {
2948: rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
2949: }
2950: for (i=0; i<nrqs; ++i) {
2951: j = pa[i];
2952: MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);
2953: }
2955: /* Send to other procs the buf size they should allocate */
2956: /* structure of sbuf2[i]/rbuf2[i]: 1 (total size to allocate) + nrow-locrow
2957: * (row sizes) */
2959: /* Receive messages*/
2960: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);
2961: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);
2962: len = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*);
2963: PetscMalloc(len,&sbuf2);
2964: req_size = (int*)(sbuf2 + nrqr);
2965: req_source = req_size + nrqr;
2966:
2967: {
2968: BSsprow **sAi = A->rows;
2969: int id,rstart = C->rmap.rstart;
2970: int *sbuf2_i,*rbuf1_i,end;
2972: for (i=0; i<nrqr; ++i) {
2973: MPI_Waitany(nrqr,r_waits1,&idx,r_status1+i);
2974: req_size[idx] = 0;
2975: rbuf1_i = rbuf1[idx];
2976: MPI_Get_count(r_status1+i,MPI_INT,&end);
2977: PetscMalloc((end+1)*sizeof(int),&sbuf2[idx]);
2978: sbuf2_i = sbuf2[idx];
2979: for (j=1; j<end; j++) {
2980: id = rbuf1_i[j] - rstart;
2981: ncols_t = (sAi[id])->length;
2982: sbuf2_i[j] = ncols_t;
2983: req_size[idx] += ncols_t;
2984: }
2985: req_source[idx] = r_status1[i].MPI_SOURCE;
2986: /* form the header */
2987: sbuf2_i[0] = req_size[idx];
2988: MPI_Isend(sbuf2_i,end,MPI_INT,req_source[idx],tag1,comm,s_waits2+i);
2989: }
2990: }
2991: PetscFree(r_status1);
2992: PetscFree(r_waits1);
2994: /* recv buffer sizes */
2995: /* Receive messages*/
2996:
2997: PetscMalloc((nrqs+1)*sizeof(int*),&rbuf3);
2998: PetscMalloc((nrqs+1)*sizeof(FLOAT*),&rbuf4);
2999: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);
3000: PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);
3001: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);
3003: for (i=0; i<nrqs; ++i) {
3004: MPI_Waitany(nrqs,r_waits2,&idx,r_status2+i);
3005: PetscMalloc((rbuf2[idx][0]+1)*sizeof(int),&rbuf3[idx]);
3006: PetscMalloc((rbuf2[idx][0]+1)*sizeof(FLOAT),&rbuf4[idx]);
3007: MPI_Irecv(rbuf3[idx],rbuf2[idx][0],MPI_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idx);
3008: MPI_Irecv(rbuf4[idx],rbuf2[idx][0],MPIU_SCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idx);
3009: }
3010: PetscFree(r_status2);
3011: PetscFree(r_waits2);
3012:
3013: /* Wait on sends1 and sends2 */
3014: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);
3015: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);
3017: if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status1);}
3018: if (nrqr) {MPI_Waitall(nrqr,s_waits2,s_status2);}
3019: PetscFree(s_status1);
3020: PetscFree(s_status2);
3021: PetscFree(s_waits1);
3022: PetscFree(s_waits2);
3024: /* Now allocate buffers for a->j, and send them off */
3025: /* structure of sbuf3[i]/rbuf3[i],sbuf4[i]/rbuf4[i]: reqsize[i] (cols resp.
3026: * vals of all req. rows; row sizes was in rbuf2; vals are of FLOAT type */
3027:
3028: PetscMalloc((nrqr+1)*sizeof(int*),&sbuf3);
3029: for (i=0,j=0; i<nrqr; i++) j += req_size[i];
3030: PetscMalloc((j+1)*sizeof(int),&sbuf3[0]);
3031: for (i=1; i<nrqr; i++) sbuf3[i] = sbuf3[i-1] + req_size[i-1];
3032:
3033: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);
3034: {
3035: int *Acol,*rbuf1_i,*sbuf3_i,rqrow,noutcols,kmax,*cols,ncols;
3036: int rstart = C->rmap.rstart;
3038: for (i=0; i<nrqr; i++) {
3039: rbuf1_i = rbuf1[i];
3040: sbuf3_i = sbuf3[i];
3041: noutcols = 0;
3042: kmax = rbuf1_i[0]; /* num. of req. rows */
3043: for (k=0,rqrow=1; k<kmax; k++,rqrow++) {
3044: Arow = A->rows[rbuf1_i[rqrow] - rstart];
3045: ncols = Arow->length;
3046: Acol = Arow->col;
3047: /* load the column indices for this row into cols*/
3048: cols = sbuf3_i + noutcols;
3049: PetscMemcpy(cols,Acol,ncols*sizeof(int));
3050: /*for (l=0; l<ncols;l++) cols[l]=Acol[l]; */ /* How is it with mappings?? */
3051: noutcols += ncols;
3052: }
3053: MPI_Isend(sbuf3_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);
3054: }
3055: }
3056: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);
3057: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);
3059: /* Allocate buffers for a->a, and send them off */
3060: /* can be optimized by conect with previous block */
3061: PetscMalloc((nrqr+1)*sizeof(FLOAT*),&sbuf4);
3062: for (i=0,j=0; i<nrqr; i++) j += req_size[i];
3063: PetscMalloc((j+1)*sizeof(FLOAT),&sbuf4[0]);
3064: for (i=1; i<nrqr; i++) sbuf4[i] = sbuf4[i-1] + req_size[i-1];
3065:
3066: PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);
3067: {
3068: FLOAT *Aval,*vals,*sbuf4_i;
3069: int rstart = C->rmap.rstart,*rbuf1_i,rqrow,noutvals,kmax,ncols;
3070:
3071:
3072: for (i=0; i<nrqr; i++) {
3073: rbuf1_i = rbuf1[i];
3074: sbuf4_i = sbuf4[i];
3075: rqrow = 1;
3076: noutvals = 0;
3077: kmax = rbuf1_i[0]; /* num of req. rows */
3078: for (k=0; k<kmax; k++,rqrow++) {
3079: Arow = A->rows[rbuf1_i[rqrow] - rstart];
3080: ncols = Arow->length;
3081: Aval = Arow->nz;
3082: /* load the column values for this row into vals*/
3083: vals = sbuf4_i+noutvals;
3084: PetscMemcpy(vals,Aval,ncols*sizeof(FLOAT));
3085: noutvals += ncols;
3086: }
3087: MPI_Isend(sbuf4_i,req_size[i],MPIU_SCALAR,req_source[i],tag3,comm,s_waits4+i);
3088: }
3089: }
3090: PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);
3091: PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);
3092: PetscFree(rbuf1);
3094: /* Form the matrix */
3096: /* create col map */
3097: len = C->cmap.N*sizeof(int)+1;
3098: PetscMalloc(len,&cmap);
3099: PetscMemzero(cmap,C->cmap.N*sizeof(int));
3100: for (j=0; j<ncol; j++) {
3101: cmap[icol[j]] = j+1;
3102: }
3103:
3104: /* Create row map / maybe I will need global rowmap but here is local rowmap*/
3105: len = C->rmap.N*sizeof(int)+1;
3106: PetscMalloc(len,&rmap);
3107: PetscMemzero(rmap,C->rmap.N*sizeof(int));
3108: for (j=0; j<nrow; j++) {
3109: rmap[irow[j]] = j;
3110: }
3112: /*
3113: Determine the number of non-zeros in the diagonal and off-diagonal
3114: portions of the matrix in order to do correct preallocation
3115: */
3117: /* first get start and end of "diagonal" columns */
3118: if (csize == PETSC_DECIDE) {
3119: nlocal = ncol/size + ((ncol % size) > rank);
3120: } else {
3121: nlocal = csize;
3122: }
3123: {
3124: int ncols,*cols,olen,dlen,thecol;
3125: int *rbuf2_i,*rbuf3_i,*sbuf1_i,row,kmax,cidx;
3126:
3127: MPI_Scan(&nlocal,&cend,1,MPI_INT,MPI_SUM,comm);
3128: cstart = cend - nlocal;
3129: if (rank == size - 1 && cend != ncol) {
3130: SETERRQ(PETSC_ERR_ARG_SIZ,"Local column sizes do not add up to total number of columns");
3131: }
3133: PetscMalloc((2*nrow+1)*sizeof(int),&d_nz);
3134: o_nz = d_nz + nrow;
3135:
3136: /* Update lens from local data */
3137: for (j=0; j<nrow; j++) {
3138: row = irow[j];
3139: proc = rtable[row];
3140: if (proc == rank) {
3141: Arow=A->rows[row-C->rmap.rstart];
3142: ncols=Arow->length;
3143: cols=Arow->col;
3144: olen=dlen=0;
3145: for (k=0; k<ncols; k++) {
3146: if ((thecol=cmap[cols[k]])) {
3147: if (cstart<thecol && thecol<=cend) dlen++; /* thecol is from 1 */
3148: else olen++;
3149: }
3150: }
3151: o_nz[j]=olen;
3152: d_nz[j]=dlen;
3153: } else d_nz[j]=o_nz[j]=0;
3154: }
3155: /* Update lens from offproc data and done waits */
3156: /* this will be much simplier after sending only appropriate columns */
3157: for (j=0; j<nrqs;j++) {
3158: MPI_Waitany(nrqs,r_waits3,&i,r_status3+j);
3159: proc = pa[i];
3160: sbuf1_i = sbuf1[proc];
3161: cidx = 0;
3162: rbuf2_i = rbuf2[i];
3163: rbuf3_i = rbuf3[i];
3164: kmax = sbuf1_i[0]; /*num of rq. rows*/
3165: for (k=1; k<=kmax; k++) {
3166: row = rmap[sbuf1_i[k]]; /* the val in the new matrix to be */
3167: for (l=0; l<rbuf2_i[k]; l++,cidx++) {
3168: if ((thecol=cmap[rbuf3_i[cidx]])) {
3169:
3170: if (cstart<thecol && thecol<=cend) d_nz[row]++; /* thecol is from 1 */
3171: else o_nz[row]++;
3172: }
3173: }
3174: }
3175: }
3176: }
3177: PetscFree(r_status3);
3178: PetscFree(r_waits3);
3179: if (nrqr) {MPI_Waitall(nrqr,s_waits3,s_status3);}
3180: PetscFree(s_status3);
3181: PetscFree(s_waits3);
3183: if (scall == MAT_INITIAL_MATRIX) {
3184: MatCreate(comm,submat);
3185: MatSetSizes(*submat,nrow,nlocal,PETSC_DECIDE,ncol);
3186: MatSetType(*submat,((PetscObject)C)->type_name);
3187: MatMPIAIJSetPreallocation(*submat,0,d_nz,0,o_nz);
3188: mat=(Mat_MPIAIJ *)((*submat)->data);
3189: matA=(Mat_SeqAIJ *)(mat->A->data);
3190: matB=(Mat_SeqAIJ *)(mat->B->data);
3191:
3192: } else {
3193: PetscTruth same;
3194: /* folowing code can be optionaly dropped for debuged versions of users
3195: * program, but I don't know PETSc option which can switch off such safety
3196: * tests - in a same way counting of o_nz,d_nz can be droped for REUSE
3197: * matrix */
3198:
3199: PetscTypeCompare((PetscObject)(*submat),MATMPIAIJ,&same);
3200: if (!same) {
3201: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong type");
3202: }
3203: if (((*submat)->rmap.n != nrow) || ((*submat)->cmap.N != ncol)) {
3204: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
3205: }
3206: mat=(Mat_MPIAIJ *)((*submat)->data);
3207: matA=(Mat_SeqAIJ *)(mat->A->data);
3208: matB=(Mat_SeqAIJ *)(mat->B->data);
3209: PetscMemcmp(matA->ilen,d_nz,nrow*sizeof(int),&same);
3210: if (!same) {
3211: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
3212: }
3213: PetscMemcmp(matB->ilen,o_nz,nrow*sizeof(int),&same);
3214: if (!same) {
3215: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
3216: }
3217: /* Initial matrix as if empty */
3218: PetscMemzero(matA->ilen,nrow*sizeof(int));
3219: PetscMemzero(matB->ilen,nrow*sizeof(int));
3220: /* Perhaps MatZeroEnteries may be better - look what it is exactly doing - I must
3221: * delete all possibly nonactual inforamtion */
3222: /*submats[i]->factor = C->factor; !!! ??? if factor will be same then I must
3223: * copy some factor information - where are thay */
3224: (*submat)->was_assembled=PETSC_FALSE;
3225: (*submat)->assembled=PETSC_FALSE;
3226:
3227: }
3228: PetscFree(d_nz);
3230: /* Assemble the matrix */
3231: /* First assemble from local rows */
3232: {
3233: int i_row,oldrow,row,ncols,*cols,*matA_j,*matB_j,ilenA,ilenB,tcol;
3234: FLOAT *vals;
3235: PetscScalar *matA_a,*matB_a;
3236:
3237: for (j=0; j<nrow; j++) {
3238: oldrow = irow[j];
3239: proc = rtable[oldrow];
3240: if (proc == rank) {
3241: row = rmap[oldrow];
3242:
3243: Arow = A->rows[oldrow-C->rmap.rstart];
3244: ncols = Arow->length;
3245: cols = Arow->col;
3246: vals = Arow->nz;
3247:
3248: i_row = matA->i[row];
3249: matA_a = matA->a + i_row;
3250: matA_j = matA->j + i_row;
3251: i_row = matB->i[row];
3252: matB_a = matB->a + i_row;
3253: matB_j = matB->j + i_row;
3254: for (k=0,ilenA=0,ilenB=0; k<ncols; k++) {
3255: if ((tcol = cmap[cols[k]])) {
3256: if (tcol<=cstart) {
3257: *matB_j++ = tcol-1;
3258: *matB_a++ = vals[k];
3259: ilenB++;
3260: } else if (tcol<=cend) {
3261: *matA_j++ = (tcol-1)-cstart;
3262: *matA_a++ = (PetscScalar)(vals[k]);
3263: ilenA++;
3264: } else {
3265: *matB_j++ = tcol-1;
3266: *matB_a++ = vals[k];
3267: ilenB++;
3268: }
3269: }
3270: }
3271: matA->ilen[row]=ilenA;
3272: matB->ilen[row]=ilenB;
3273:
3274: }
3275: }
3276: }
3278: /* Now assemble the off proc rows*/
3279: {
3280: int *sbuf1_i,*rbuf2_i,*rbuf3_i,cidx,kmax,row,i_row;
3281: int *matA_j,*matB_j,lmax,tcol,ilenA,ilenB;
3282: PetscScalar *matA_a,*matB_a;
3283: FLOAT *rbuf4_i;
3285: for (j=0; j<nrqs; j++) {
3286: MPI_Waitany(nrqs,r_waits4,&i,r_status4+j);
3287: proc = pa[i];
3288: sbuf1_i = sbuf1[proc];
3289:
3290: cidx = 0;
3291: rbuf2_i = rbuf2[i];
3292: rbuf3_i = rbuf3[i];
3293: rbuf4_i = rbuf4[i];
3294: kmax = sbuf1_i[0];
3295: for (k=1; k<=kmax; k++) {
3296: row = rmap[sbuf1_i[k]];
3297:
3298: i_row = matA->i[row];
3299: matA_a = matA->a + i_row;
3300: matA_j = matA->j + i_row;
3301: i_row = matB->i[row];
3302: matB_a = matB->a + i_row;
3303: matB_j = matB->j + i_row;
3304:
3305: lmax = rbuf2_i[k];
3306: for (l=0,ilenA=0,ilenB=0; l<lmax; l++,cidx++) {
3307: if ((tcol = cmap[rbuf3_i[cidx]])) {
3308: if (tcol<=cstart) {
3309: *matB_j++ = tcol-1;
3310: *matB_a++ = (PetscScalar)(rbuf4_i[cidx]);;
3311: ilenB++;
3312: } else if (tcol<=cend) {
3313: *matA_j++ = (tcol-1)-cstart;
3314: *matA_a++ = (PetscScalar)(rbuf4_i[cidx]);
3315: ilenA++;
3316: } else {
3317: *matB_j++ = tcol-1;
3318: *matB_a++ = (PetscScalar)(rbuf4_i[cidx]);
3319: ilenB++;
3320: }
3321: }
3322: }
3323: matA->ilen[row]=ilenA;
3324: matB->ilen[row]=ilenB;
3325: }
3326: }
3327: }
3329: PetscFree(r_status4);
3330: PetscFree(r_waits4);
3331: if (nrqr) {MPI_Waitall(nrqr,s_waits4,s_status4);}
3332: PetscFree(s_waits4);
3333: PetscFree(s_status4);
3335: /* Restore the indices */
3336: ISRestoreIndices(isrow,&irow);
3337: ISRestoreIndices(iscol,&icol);
3339: /* Destroy allocated memory */
3340: PetscFree(rtable);
3341: PetscFree(w1);
3342: PetscFree(pa);
3344: PetscFree(sbuf1);
3345: PetscFree(rbuf2[0]);
3346: PetscFree(rbuf2);
3347: for (i=0; i<nrqr; ++i) {
3348: PetscFree(sbuf2[i]);
3349: }
3350: for (i=0; i<nrqs; ++i) {
3351: PetscFree(rbuf3[i]);
3352: PetscFree(rbuf4[i]);
3353: }
3355: PetscFree(sbuf2);
3356: PetscFree(rbuf3);
3357: PetscFree(rbuf4);
3358: PetscFree(sbuf3[0]);
3359: PetscFree(sbuf3);
3360: PetscFree(sbuf4[0]);
3361: PetscFree(sbuf4);
3362:
3363: PetscFree(cmap);
3364: PetscFree(rmap);
3367: MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);
3368: MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);
3371: return(0);
3372: }