Actual source code: mpibdiag.c
1: #define PETSCMAT_DLL
3: /*
4: The basic matrix operations for the Block diagonal parallel
5: matrices.
6: */
7: #include src/mat/impls/bdiag/mpi/mpibdiag.h
11: PetscErrorCode MatSetValues_MPIBDiag(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
12: {
13: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
15: PetscInt i,j,row,rstart = mat->rmap.rstart,rend = mat->rmap.rend;
16: PetscTruth roworiented = mbd->roworiented;
19: for (i=0; i<m; i++) {
20: if (idxm[i] < 0) continue;
21: if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
22: if (idxm[i] >= rstart && idxm[i] < rend) {
23: row = idxm[i] - rstart;
24: for (j=0; j<n; j++) {
25: if (idxn[j] < 0) continue;
26: if (idxn[j] >= mat->cmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
27: if (roworiented) {
28: MatSetValues(mbd->A,1,&row,1,&idxn[j],v+i*n+j,addv);
29: } else {
30: MatSetValues(mbd->A,1,&row,1,&idxn[j],v+i+j*m,addv);
31: }
32: }
33: } else {
34: if (!mbd->donotstash) {
35: if (roworiented) {
36: MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);
37: } else {
38: MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);
39: }
40: }
41: }
42: }
43: return(0);
44: }
48: PetscErrorCode MatGetValues_MPIBDiag(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
49: {
50: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
52: PetscInt i,j,row,rstart = mat->rmap.rstart,rend = mat->rmap.rend;
55: for (i=0; i<m; i++) {
56: if (idxm[i] < 0) continue; /* SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
57: if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
58: if (idxm[i] >= rstart && idxm[i] < rend) {
59: row = idxm[i] - rstart;
60: for (j=0; j<n; j++) {
61: if (idxn[j] < 0) continue; /* SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
62: if (idxn[j] >= mat->cmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
63: MatGetValues(mbd->A,1,&row,1,&idxn[j],v+i*n+j);
64: }
65: } else {
66: SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
67: }
68: }
69: return(0);
70: }
74: PetscErrorCode MatAssemblyBegin_MPIBDiag(Mat mat,MatAssemblyType mode)
75: {
76: MPI_Comm comm = ((PetscObject)mat)->comm;
78: PetscInt nstash,reallocs;
79: InsertMode addv;
82: MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
83: if (addv == (ADD_VALUES|INSERT_VALUES)) {
84: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
85: }
86: mat->insertmode = addv; /* in case this processor had no cache */
87: MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap.range);
88: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
89: PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
90: return(0);
91: }
95: PetscErrorCode MatAssemblyEnd_MPIBDiag(Mat mat,MatAssemblyType mode)
96: {
97: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
98: Mat_SeqBDiag *mlocal;
100: PetscMPIInt n;
101: PetscInt i,*row,*col;
102: PetscInt *tmp1,*tmp2,len,ict,Mblock,Nblock,flg,j,rstart,ncols;
103: PetscScalar *val;
104: InsertMode addv = mat->insertmode;
108: while (1) {
109: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
110: if (!flg) break;
111:
112: for (i=0; i<n;) {
113: /* Now identify the consecutive vals belonging to the same row */
114: for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
115: if (j < n) ncols = j-i;
116: else ncols = n-i;
117: /* Now assemble all these values with a single function call */
118: MatSetValues_MPIBDiag(mat,1,row+i,ncols,col+i,val+i,addv);
119: i = j;
120: }
121: }
122: MatStashScatterEnd_Private(&mat->stash);
124: MatAssemblyBegin(mbd->A,mode);
125: MatAssemblyEnd(mbd->A,mode);
127: /* Fix main diagonal location and determine global diagonals */
128: mlocal = (Mat_SeqBDiag*)mbd->A->data;
129: Mblock = mat->rmap.N/mat->rmap.bs; Nblock = mat->cmap.N/mat->rmap.bs;
130: len = Mblock + Nblock + 1; /* add 1 to prevent 0 malloc */
131: PetscMalloc(2*len*sizeof(PetscInt),&tmp1);
132: tmp2 = tmp1 + len;
133: PetscMemzero(tmp1,2*len*sizeof(PetscInt));
134: mlocal->mainbd = -1;
135: for (i=0; i<mlocal->nd; i++) {
136: if (mlocal->diag[i] + mbd->brstart == 0) mlocal->mainbd = i;
137: tmp1[mlocal->diag[i] + mbd->brstart + Mblock] = 1;
138: }
139: MPI_Allreduce(tmp1,tmp2,len,MPIU_INT,MPI_SUM,((PetscObject)mat)->comm);
140: ict = 0;
141: for (i=0; i<len; i++) {
142: if (tmp2[i]) {
143: mbd->gdiag[ict] = i - Mblock;
144: ict++;
145: }
146: }
147: mbd->gnd = ict;
148: PetscFree(tmp1);
150: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
151: MatSetUpMultiply_MPIBDiag(mat);
152: }
153: return(0);
154: }
158: PetscErrorCode MatZeroEntries_MPIBDiag(Mat A)
159: {
160: Mat_MPIBDiag *l = (Mat_MPIBDiag*)A->data;
164: MatZeroEntries(l->A);
165: return(0);
166: }
168: /* again this uses the same basic stratagy as in the assembly and
169: scatter create routines, we should try to do it systematically
170: if we can figure out the proper level of generality. */
172: /* the code does not do the diagonal entries correctly unless the
173: matrix is square and the column and row owerships are identical.
174: This is a BUG. The only way to fix it seems to be to access
175: aij->A and aij->B directly and not through the MatZeroRows()
176: routine.
177: */
181: PetscErrorCode MatZeroRows_MPIBDiag(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
182: {
183: Mat_MPIBDiag *l = (Mat_MPIBDiag*)A->data;
185: PetscMPIInt n,imdex,size = l->size,rank = l->rank,tag = ((PetscObject)A)->tag;
186: PetscInt i,*owners = A->rmap.range;
187: PetscInt *nprocs,j,idx,nsends;
188: PetscInt nmax,*svalues,*starts,*owner,nrecvs;
189: PetscInt *rvalues,count,base,slen,*source;
190: PetscInt *lens,*lrows,*values;
191: MPI_Comm comm = ((PetscObject)A)->comm;
192: MPI_Request *send_waits,*recv_waits;
193: MPI_Status recv_status,*send_status;
194: PetscTruth found;
197: /* first count number of contributors to each processor */
198: PetscMalloc(2*size*sizeof(PetscInt),&nprocs);
199: PetscMemzero(nprocs,2*size*sizeof(PetscInt));
200: PetscMalloc((N+1)*sizeof(PetscInt),&owner); /* see note*/
201: for (i=0; i<N; i++) {
202: idx = rows[i];
203: found = PETSC_FALSE;
204: for (j=0; j<size; j++) {
205: if (idx >= owners[j] && idx < owners[j+1]) {
206: nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
207: }
208: }
209: if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
210: }
211: nsends = 0; for (i=0; i<size; i++) {nsends += nprocs[2*i+1];}
213: /* inform other processors of number of messages and max length*/
214: PetscMaxSum(comm,nprocs,&nmax,&nrecvs);
216: /* post receives: */
217: PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);
218: PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);
219: for (i=0; i<nrecvs; i++) {
220: MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
221: }
223: /* do sends:
224: 1) starts[i] gives the starting index in svalues for stuff going to
225: the ith processor
226: */
227: PetscMalloc((N+1)*sizeof(PetscInt),&svalues);
228: PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);
229: PetscMalloc((size+1)*sizeof(PetscInt),&starts);
230: starts[0] = 0;
231: for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
232: for (i=0; i<N; i++) {
233: svalues[starts[owner[i]]++] = rows[i];
234: }
236: starts[0] = 0;
237: for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
238: count = 0;
239: for (i=0; i<size; i++) {
240: if (nprocs[2*i+1]) {
241: MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);
242: }
243: }
244: PetscFree(starts);
246: base = owners[rank];
248: /* wait on receives */
249: PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);
250: source = lens + nrecvs;
251: count = nrecvs;
252: slen = 0;
253: while (count) {
254: MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
255: /* unpack receives into our local space */
256: MPI_Get_count(&recv_status,MPIU_INT,&n);
257: source[imdex] = recv_status.MPI_SOURCE;
258: lens[imdex] = n;
259: slen += n;
260: count--;
261: }
262: PetscFree(recv_waits);
263:
264: /* move the data into the send scatter */
265: PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);
266: count = 0;
267: for (i=0; i<nrecvs; i++) {
268: values = rvalues + i*nmax;
269: for (j=0; j<lens[i]; j++) {
270: lrows[count++] = values[j] - base;
271: }
272: }
273: PetscFree(rvalues);
274: PetscFree(lens);
275: PetscFree(owner);
276: PetscFree(nprocs);
277:
278: /* actually zap the local rows */
279: MatZeroRows(l->A,slen,lrows,diag);
280: PetscFree(lrows);
282: /* wait on sends */
283: if (nsends) {
284: PetscMalloc(nsends*sizeof(MPI_Status),&send_status);
285: MPI_Waitall(nsends,send_waits,send_status);
286: PetscFree(send_status);
287: }
288: PetscFree(send_waits);
289: PetscFree(svalues);
291: return(0);
292: }
296: PetscErrorCode MatMult_MPIBDiag(Mat mat,Vec xx,Vec yy)
297: {
298: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
302: VecScatterBegin(mbd->Mvctx,xx,mbd->lvec,INSERT_VALUES,SCATTER_FORWARD);
303: VecScatterEnd(mbd->Mvctx,xx,mbd->lvec,INSERT_VALUES,SCATTER_FORWARD);
304: (*mbd->A->ops->mult)(mbd->A,mbd->lvec,yy);
305: return(0);
306: }
310: PetscErrorCode MatMultAdd_MPIBDiag(Mat mat,Vec xx,Vec yy,Vec zz)
311: {
312: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
316: VecScatterBegin(mbd->Mvctx,xx,mbd->lvec,INSERT_VALUES,SCATTER_FORWARD);
317: VecScatterEnd(mbd->Mvctx,xx,mbd->lvec,INSERT_VALUES,SCATTER_FORWARD);
318: (*mbd->A->ops->multadd)(mbd->A,mbd->lvec,yy,zz);
319: return(0);
320: }
324: PetscErrorCode MatMultTranspose_MPIBDiag(Mat A,Vec xx,Vec yy)
325: {
326: Mat_MPIBDiag *a = (Mat_MPIBDiag*)A->data;
328: PetscScalar zero = 0.0;
331: VecSet(yy,zero);
332: (*a->A->ops->multtranspose)(a->A,xx,a->lvec);
333: VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
334: VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
335: return(0);
336: }
340: PetscErrorCode MatMultTransposeAdd_MPIBDiag(Mat A,Vec xx,Vec yy,Vec zz)
341: {
342: Mat_MPIBDiag *a = (Mat_MPIBDiag*)A->data;
346: VecCopy(yy,zz);
347: (*a->A->ops->multtranspose)(a->A,xx,a->lvec);
348: VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
349: VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
350: return(0);
351: }
355: PetscErrorCode MatGetInfo_MPIBDiag(Mat matin,MatInfoType flag,MatInfo *info)
356: {
357: Mat_MPIBDiag *mat = (Mat_MPIBDiag*)matin->data;
359: PetscReal isend[5],irecv[5];
362: info->block_size = (PetscReal)mat->A->rmap.bs;
363: MatGetInfo(mat->A,MAT_LOCAL,info);
364: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
365: isend[3] = info->memory; isend[4] = info->mallocs;
366: if (flag == MAT_LOCAL) {
367: info->nz_used = isend[0];
368: info->nz_allocated = isend[1];
369: info->nz_unneeded = isend[2];
370: info->memory = isend[3];
371: info->mallocs = isend[4];
372: } else if (flag == MAT_GLOBAL_MAX) {
373: MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,((PetscObject)matin)->comm);
374: info->nz_used = irecv[0];
375: info->nz_allocated = irecv[1];
376: info->nz_unneeded = irecv[2];
377: info->memory = irecv[3];
378: info->mallocs = irecv[4];
379: } else if (flag == MAT_GLOBAL_SUM) {
380: MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,((PetscObject)matin)->comm);
381: info->nz_used = irecv[0];
382: info->nz_allocated = irecv[1];
383: info->nz_unneeded = irecv[2];
384: info->memory = irecv[3];
385: info->mallocs = irecv[4];
386: }
387: info->rows_global = (double)matin->rmap.N;
388: info->columns_global = (double)matin->cmap.N;
389: info->rows_local = (double)matin->rmap.n;
390: info->columns_local = (double)matin->cmap.N;
391: return(0);
392: }
396: PetscErrorCode MatGetDiagonal_MPIBDiag(Mat mat,Vec v)
397: {
399: Mat_MPIBDiag *A = (Mat_MPIBDiag*)mat->data;
402: MatGetDiagonal(A->A,v);
403: return(0);
404: }
408: PetscErrorCode MatDestroy_MPIBDiag(Mat mat)
409: {
410: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
412: #if defined(PETSC_USE_LOG)
413: Mat_SeqBDiag *ms = (Mat_SeqBDiag*)mbd->A->data;
416: PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, BSize=%D, NDiag=%D",mat->rmap.N,mat->cmap.N,mat->rmap.bs,ms->nd);
417: #else
419: #endif
420: MatStashDestroy_Private(&mat->stash);
421: PetscFree(mbd->gdiag);
422: MatDestroy(mbd->A);
423: if (mbd->lvec) {VecDestroy(mbd->lvec);}
424: if (mbd->Mvctx) {VecScatterDestroy(mbd->Mvctx);}
425: PetscFree(mbd);
427: PetscObjectChangeTypeName((PetscObject)mat,0);
428: PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);
429: PetscObjectComposeFunction((PetscObject)mat,"MatMPIBDiagSetPreallocation_C","",PETSC_NULL);
430: return(0);
431: }
436: static PetscErrorCode MatView_MPIBDiag_Binary(Mat mat,PetscViewer viewer)
437: {
438: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
442: if (mbd->size == 1) {
443: MatView(mbd->A,viewer);
444: } else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
445: return(0);
446: }
450: static PetscErrorCode MatView_MPIBDiag_ASCIIorDraw(Mat mat,PetscViewer viewer)
451: {
452: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
453: PetscErrorCode ierr;
454: PetscMPIInt size = mbd->size,rank = mbd->rank;
455: PetscInt i;
456: PetscTruth iascii,isdraw;
457: PetscViewer sviewer;
458: PetscViewerFormat format;
461: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
462: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
463: if (iascii) {
464: PetscViewerGetFormat(viewer,&format);
465: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
466: PetscInt nline = PetscMin(10,mbd->gnd),k,nk,np;
467: PetscViewerASCIIPrintf(viewer," block size=%D, total number of diagonals=%D\n",mat->rmap.bs,mbd->gnd);
468: nk = (mbd->gnd-1)/nline + 1;
469: for (k=0; k<nk; k++) {
470: PetscViewerASCIIPrintf(viewer," global diag numbers:");
471: np = PetscMin(nline,mbd->gnd - nline*k);
472: for (i=0; i<np; i++) {
473: PetscViewerASCIIPrintf(viewer," %D",mbd->gdiag[i+nline*k]);
474: }
475: PetscViewerASCIIPrintf(viewer,"\n");
476: }
477: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
478: MatInfo info;
479: MPI_Comm_rank(((PetscObject)mat)->comm,&rank);
480: MatGetInfo(mat,MAT_LOCAL,&info);
481: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap.N,
482: (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);
483: PetscViewerFlush(viewer);
484: VecScatterView(mbd->Mvctx,viewer);
485: }
486: return(0);
487: }
488: }
490: if (isdraw) {
491: PetscDraw draw;
492: PetscTruth isnull;
493: PetscViewerDrawGetDraw(viewer,0,&draw);
494: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
495: }
497: if (size == 1) {
498: MatView(mbd->A,viewer);
499: } else {
500: /* assemble the entire matrix onto first processor. */
501: Mat A;
502: PetscInt M = mat->rmap.N,N = mat->cmap.N,m,row,nz,*cols;
503: PetscScalar *vals;
505: /* Here we are constructing a temporary matrix, so we will explicitly set the type to MPIBDiag */
506: MatCreate(((PetscObject)mat)->comm,&A);
507: if (!rank) {
508: MatSetSizes(A,M,N,M,N);
509: MatSetType(A,MATMPIBDIAG);
510: MatMPIBDiagSetPreallocation(A,mbd->gnd,mbd->A->rmap.bs,mbd->gdiag,PETSC_NULL);
511: } else {
512: MatSetSizes(A,0,0,M,N);
513: MatSetType(A,MATMPIBDIAG);
514: MatMPIBDiagSetPreallocation(A,0,mbd->A->rmap.bs,PETSC_NULL,PETSC_NULL);
515: }
516: PetscLogObjectParent(mat,A);
518: /* Copy the matrix ... This isn't the most efficient means,
519: but it's quick for now */
520: row = mat->rmap.rstart;
521: m = mbd->A->rmap.N;
522: for (i=0; i<m; i++) {
523: MatGetRow_MPIBDiag(mat,row,&nz,&cols,&vals);
524: MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);
525: MatRestoreRow_MPIBDiag(mat,row,&nz,&cols,&vals);
526: row++;
527: }
528: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
529: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
530: PetscViewerGetSingleton(viewer,&sviewer);
531: if (!rank) {
532: MatView(((Mat_MPIBDiag*)(A->data))->A,sviewer);
533: }
534: PetscViewerRestoreSingleton(viewer,&sviewer);
535: PetscViewerFlush(viewer);
536: MatDestroy(A);
537: }
538: return(0);
539: }
543: PetscErrorCode MatView_MPIBDiag(Mat mat,PetscViewer viewer)
544: {
546: PetscTruth iascii,isdraw,isbinary;
549: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
550: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
551: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
552: if (iascii || isdraw) {
553: MatView_MPIBDiag_ASCIIorDraw(mat,viewer);
554: } else if (isbinary) {
555: MatView_MPIBDiag_Binary(mat,viewer);
556: } else {
557: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPIBdiag matrices",((PetscObject)viewer)->type_name);
558: }
559: return(0);
560: }
564: PetscErrorCode MatSetOption_MPIBDiag(Mat A,MatOption op,PetscTruth flg)
565: {
566: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)A->data;
569: switch (op) {
570: case MAT_NEW_NONZERO_LOCATIONS:
571: case MAT_NEW_NONZERO_LOCATION_ERR:
572: case MAT_NEW_NONZERO_ALLOCATION_ERR:
573: case MAT_NEW_DIAGONALS:
574: MatSetOption(mbd->A,op,flg);
575: break;
576: case MAT_ROW_ORIENTED:
577: mbd->roworiented = flg;
578: MatSetOption(mbd->A,op,flg);
579: break;
580: case MAT_IGNORE_OFF_PROC_ENTRIES:
581: mbd->donotstash = flg;
582: break;
583: case MAT_SYMMETRIC:
584: case MAT_STRUCTURALLY_SYMMETRIC:
585: case MAT_HERMITIAN:
586: case MAT_SYMMETRY_ETERNAL:
587: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
588: break;
589: default:
590: SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
591: }
592: return(0);
593: }
598: PetscErrorCode MatGetRow_MPIBDiag(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
599: {
600: Mat_MPIBDiag *mat = (Mat_MPIBDiag*)matin->data;
602: PetscInt lrow;
605: if (row < matin->rmap.rstart || row >= matin->rmap.rend) SETERRQ(PETSC_ERR_SUP,"only for local rows")
606: lrow = row - matin->rmap.rstart;
607: MatGetRow_SeqBDiag(mat->A,lrow,nz,idx,v);
608: return(0);
609: }
613: PetscErrorCode MatRestoreRow_MPIBDiag(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
614: {
615: Mat_MPIBDiag *mat = (Mat_MPIBDiag*)matin->data;
617: PetscInt lrow;
620: lrow = row - matin->rmap.rstart;
621: MatRestoreRow_SeqBDiag(mat->A,lrow,nz,idx,v);
622: return(0);
623: }
628: PetscErrorCode MatNorm_MPIBDiag(Mat A,NormType type,PetscReal *nrm)
629: {
630: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)A->data;
631: Mat_SeqBDiag *a = (Mat_SeqBDiag*)mbd->A->data;
632: PetscReal sum = 0.0;
634: PetscInt d,i,nd = a->nd,bs = A->rmap.bs,len;
635: PetscScalar *dv;
638: if (type == NORM_FROBENIUS) {
639: for (d=0; d<nd; d++) {
640: dv = a->diagv[d];
641: len = a->bdlen[d]*bs*bs;
642: for (i=0; i<len; i++) {
643: #if defined(PETSC_USE_COMPLEX)
644: sum += PetscRealPart(PetscConj(dv[i])*dv[i]);
645: #else
646: sum += dv[i]*dv[i];
647: #endif
648: }
649: }
650: MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);
651: *nrm = sqrt(*nrm);
652: PetscLogFlops(2*A->rmap.n*A->rmap.N);
653: } else if (type == NORM_1) { /* max column norm */
654: PetscReal *tmp,*tmp2;
655: PetscInt j;
656: PetscMalloc((mbd->A->cmap.n+1)*sizeof(PetscReal),&tmp);
657: PetscMalloc((mbd->A->cmap.n+1)*sizeof(PetscReal),&tmp2);
658: MatNorm_SeqBDiag_Columns(mbd->A,tmp,mbd->A->cmap.n);
659: *nrm = 0.0;
660: MPI_Allreduce(tmp,tmp2,mbd->A->cmap.n,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);
661: for (j=0; j<mbd->A->cmap.n; j++) {
662: if (tmp2[j] > *nrm) *nrm = tmp2[j];
663: }
664: PetscFree(tmp);
665: PetscFree(tmp2);
666: } else if (type == NORM_INFINITY) { /* max row norm */
667: PetscReal normtemp;
668: MatNorm(mbd->A,type,&normtemp);
669: MPI_Allreduce(&normtemp,nrm,1,MPIU_REAL,MPI_MAX,((PetscObject)A)->comm);
670: }
671: return(0);
672: }
676: PetscErrorCode MatScale_MPIBDiag(Mat A,PetscScalar alpha)
677: {
679: Mat_MPIBDiag *a = (Mat_MPIBDiag*)A->data;
682: MatScale_SeqBDiag(a->A,alpha);
683: return(0);
684: }
688: PetscErrorCode MatSetUpPreallocation_MPIBDiag(Mat A)
689: {
693: MatMPIBDiagSetPreallocation(A,PETSC_DEFAULT,PETSC_DEFAULT,0,0);
694: return(0);
695: }
697: /* -------------------------------------------------------------------*/
699: static struct _MatOps MatOps_Values = {MatSetValues_MPIBDiag,
700: MatGetRow_MPIBDiag,
701: MatRestoreRow_MPIBDiag,
702: MatMult_MPIBDiag,
703: /* 4*/ MatMultAdd_MPIBDiag,
704: MatMultTranspose_MPIBDiag,
705: MatMultTransposeAdd_MPIBDiag,
706: 0,
707: 0,
708: 0,
709: /*10*/ 0,
710: 0,
711: 0,
712: 0,
713: 0,
714: /*15*/ MatGetInfo_MPIBDiag,
715: 0,
716: MatGetDiagonal_MPIBDiag,
717: 0,
718: MatNorm_MPIBDiag,
719: /*20*/ MatAssemblyBegin_MPIBDiag,
720: MatAssemblyEnd_MPIBDiag,
721: 0,
722: MatSetOption_MPIBDiag,
723: MatZeroEntries_MPIBDiag,
724: /*25*/ MatZeroRows_MPIBDiag,
725: 0,
726: 0,
727: 0,
728: 0,
729: /*30*/ MatSetUpPreallocation_MPIBDiag,
730: 0,
731: 0,
732: 0,
733: 0,
734: /*35*/ 0,
735: 0,
736: 0,
737: 0,
738: 0,
739: /*40*/ 0,
740: 0,
741: 0,
742: MatGetValues_MPIBDiag,
743: 0,
744: /*45*/ 0,
745: MatScale_MPIBDiag,
746: 0,
747: 0,
748: 0,
749: /*50*/ 0,
750: 0,
751: 0,
752: 0,
753: 0,
754: /*55*/ 0,
755: 0,
756: 0,
757: 0,
758: 0,
759: /*60*/ 0,
760: MatDestroy_MPIBDiag,
761: MatView_MPIBDiag,
762: 0,
763: 0,
764: /*65*/ 0,
765: 0,
766: 0,
767: 0,
768: 0,
769: /*70*/ 0,
770: 0,
771: 0,
772: 0,
773: 0,
774: /*75*/ 0,
775: 0,
776: 0,
777: 0,
778: 0,
779: /*80*/ 0,
780: 0,
781: 0,
782: 0,
783: MatLoad_MPIBDiag,
784: /*85*/ 0,
785: 0,
786: 0,
787: 0,
788: 0,
789: /*90*/ 0,
790: 0,
791: 0,
792: 0,
793: 0,
794: /*95*/ 0,
795: 0,
796: 0,
797: 0};
802: PetscErrorCode MatGetDiagonalBlock_MPIBDiag(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *a)
803: {
804: Mat_MPIBDiag *matin = (Mat_MPIBDiag *)A->data;
806: PetscInt lrows,lcols,rstart,rend;
807: IS localc,localr;
810: MatGetLocalSize(A,&lrows,&lcols);
811: MatGetOwnershipRange(A,&rstart,&rend);
812: ISCreateStride(PETSC_COMM_SELF,lrows,rstart,1,&localc);
813: ISCreateStride(PETSC_COMM_SELF,lrows,0,1,&localr);
814: MatGetSubMatrix(matin->A,localr,localc,PETSC_DECIDE,reuse,a);
815: ISDestroy(localr);
816: ISDestroy(localc);
818: *iscopy = PETSC_TRUE;
819: return(0);
820: }
826: PetscErrorCode MatMPIBDiagSetPreallocation_MPIBDiag(Mat B,PetscInt nd,PetscInt bs,PetscInt *diag,PetscScalar **diagv)
827: {
828: Mat_MPIBDiag *b;
830: PetscInt i,k,*ldiag,len,nd2;
831: PetscScalar **ldiagv = 0;
832: PetscTruth flg2;
835: B->preallocated = PETSC_TRUE;
836: if (bs == PETSC_DEFAULT) bs = 1;
837: if (nd == PETSC_DEFAULT) nd = 0;
838: PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
839: PetscOptionsGetInt(PETSC_NULL,"-mat_bdiag_ndiag",&nd,PETSC_NULL);
840: PetscOptionsHasName(PETSC_NULL,"-mat_bdiag_diags",&flg2);
841: if (nd && !diag) {
842: PetscMalloc(nd*sizeof(PetscInt),&diag);
843: nd2 = nd;
844: PetscOptionsGetIntArray(PETSC_NULL,"-mat_bdiag_dvals",diag,&nd2,PETSC_NULL);
845: if (nd2 != nd) {
846: SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible number of diags and diagonal vals");
847: }
848: } else if (flg2) {
849: SETERRQ(PETSC_ERR_ARG_WRONG,"Must specify number of diagonals with -mat_bdiag_ndiag");
850: }
852: if (bs <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive");
854: B->rmap.bs = B->cmap.bs = bs;
856: PetscMapSetUp(&B->rmap);
857: PetscMapSetUp(&B->cmap);
859: if ((B->cmap.N%bs)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size - bad column number");
860: if ((B->rmap.N%bs)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size - bad local row number");
861: if ((B->rmap.N%bs)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size - bad global row number");
864: b = (Mat_MPIBDiag*)B->data;
865: b->gnd = nd;
866: b->brstart = (B->rmap.rstart)/bs;
867: b->brend = (B->rmap.rend)/bs;
870: /* Determine local diagonals; for now, assume global rows = global cols */
871: /* These are sorted in MatCreateSeqBDiag */
872: PetscMalloc((nd+1)*sizeof(PetscInt),&ldiag);
873: len = B->rmap.N/bs + B->cmap.N/bs + 1;
874: PetscMalloc(len*sizeof(PetscInt),&b->gdiag);
875: k = 0;
876: if (diagv) {
877: PetscMalloc((nd+1)*sizeof(PetscScalar*),&ldiagv);
878: }
879: for (i=0; i<nd; i++) {
880: b->gdiag[i] = diag[i];
881: if (diag[i] > 0) { /* lower triangular */
882: if (diag[i] < b->brend) {
883: ldiag[k] = diag[i] - b->brstart;
884: if (diagv) ldiagv[k] = diagv[i];
885: k++;
886: }
887: } else { /* upper triangular */
888: if (B->rmap.N/bs - diag[i] > B->cmap.N/bs) {
889: if (B->rmap.N/bs + diag[i] > b->brstart) {
890: ldiag[k] = diag[i] - b->brstart;
891: if (diagv) ldiagv[k] = diagv[i];
892: k++;
893: }
894: } else {
895: if (B->rmap.N/bs > b->brstart) {
896: ldiag[k] = diag[i] - b->brstart;
897: if (diagv) ldiagv[k] = diagv[i];
898: k++;
899: }
900: }
901: }
902: }
904: /* Form local matrix */
905: MatCreate(PETSC_COMM_SELF,&b->A);
906: MatSetSizes(b->A,B->rmap.n,B->cmap.N,B->rmap.n,B->cmap.N);
907: MatSetType(b->A,MATSEQBDIAG);
908: MatSeqBDiagSetPreallocation(b->A,k,bs,ldiag,ldiagv);
909: PetscLogObjectParent(B,b->A);
910: PetscFree(ldiag);
911: PetscFree(ldiagv);
913: return(0);
914: }
917: /*MC
918: MATMPIBDIAG - MATMPIBDIAG = "mpibdiag" - A matrix type to be used for distributed block diagonal matrices.
920: Options Database Keys:
921: . -mat_type mpibdiag - sets the matrix type to "mpibdiag" during a call to MatSetFromOptions()
923: Level: beginner
925: .seealso: MatCreateMPIBDiag
926: M*/
931: PetscErrorCode MatCreate_MPIBDiag(Mat B)
932: {
933: Mat_MPIBDiag *b;
937: PetscNewLog(B,Mat_MPIBDiag,&b);
938: B->data = (void*)b;
939: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
940: B->factor = 0;
941: B->mapping = 0;
943: B->insertmode = NOT_SET_VALUES;
944: MPI_Comm_rank(((PetscObject)B)->comm,&b->rank);
945: MPI_Comm_size(((PetscObject)B)->comm,&b->size);
947: /* build cache for off array entries formed */
948: MatStashCreate_Private(((PetscObject)B)->comm,1,&B->stash);
949: b->donotstash = PETSC_FALSE;
951: /* stuff used for matrix-vector multiply */
952: b->lvec = 0;
953: b->Mvctx = 0;
955: /* used for MatSetValues() input */
956: b->roworiented = PETSC_TRUE;
958: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetDiagonalBlock_C",
959: "MatGetDiagonalBlock_MPIBDiag",
960: MatGetDiagonalBlock_MPIBDiag);
961: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIBDiagSetPreallocation_C",
962: "MatMPIBDiagSetPreallocation_MPIBDiag",
963: MatMPIBDiagSetPreallocation_MPIBDiag);
964: PetscObjectChangeTypeName((PetscObject)B,MATMPIBDIAG);
965: return(0);
966: }
969: /*MC
970: MATBDIAG - MATBDIAG = "bdiag" - A matrix type to be used for block diagonal matrices.
972: This matrix type is identical to MATSEQBDIAG when constructed with a single process communicator,
973: and MATMPIBDIAG otherwise.
975: Options Database Keys:
976: . -mat_type bdiag - sets the matrix type to "bdiag" during a call to MatSetFromOptions()
978: Level: beginner
980: .seealso: MatCreateMPIBDiag,MATSEQBDIAG,MATMPIBDIAG
981: M*/
986: PetscErrorCode MatCreate_BDiag(Mat A)
987: {
989: PetscMPIInt size;
992: MPI_Comm_size(((PetscObject)A)->comm,&size);
993: if (size == 1) {
994: MatSetType(A,MATSEQBDIAG);
995: } else {
996: MatSetType(A,MATMPIBDIAG);
997: }
998: return(0);
999: }
1004: /*@C
1005: MatMPIBDiagSetPreallocation -
1007: Collective on Mat
1009: Input Parameters:
1010: + A - the matrix
1011: . nd - number of block diagonals (global) (optional)
1012: . bs - each element of a diagonal is an bs x bs dense matrix
1013: . diag - optional array of block diagonal numbers (length nd).
1014: For a matrix element A[i,j], where i=row and j=column, the
1015: diagonal number is
1016: $ diag = i/bs - j/bs (integer division)
1017: Set diag=PETSC_NULL on input for PETSc to dynamically allocate memory as
1018: needed (expensive).
1019: - diagv - pointer to actual diagonals (in same order as diag array),
1020: if allocated by user. Otherwise, set diagv=PETSC_NULL on input for PETSc
1021: to control memory allocation.
1024: Options Database Keys:
1025: . -mat_block_size <bs> - Sets blocksize
1026: . -mat_bdiag_diags <s1,s2,s3,...> - Sets diagonal numbers
1028: Notes:
1029: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
1030: than it must be used on all processors that share the object for that argument.
1032: The parallel matrix is partitioned across the processors by rows, where
1033: each local rectangular matrix is stored in the uniprocessor block
1034: diagonal format. See the users manual for further details.
1036: The user MUST specify either the local or global numbers of rows
1037: (possibly both).
1039: The case bs=1 (conventional diagonal storage) is implemented as
1040: a special case.
1042: Fortran Notes:
1043: Fortran programmers cannot set diagv; this variable is ignored.
1045: Level: intermediate
1047: .keywords: matrix, block, diagonal, parallel, sparse
1049: .seealso: MatCreate(), MatCreateSeqBDiag(), MatSetValues()
1050: @*/
1051: PetscErrorCode MatMPIBDiagSetPreallocation(Mat B,PetscInt nd,PetscInt bs,const PetscInt diag[],PetscScalar *diagv[])
1052: {
1053: PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
1056: PetscObjectQueryFunction((PetscObject)B,"MatMPIBDiagSetPreallocation_C",(void (**)(void))&f);
1057: if (f) {
1058: (*f)(B,nd,bs,diag,diagv);
1059: }
1060: return(0);
1061: }
1065: /*@C
1066: MatCreateMPIBDiag - Creates a sparse parallel matrix in MPIBDiag format.
1068: Collective on MPI_Comm
1070: Input Parameters:
1071: + comm - MPI communicator
1072: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1073: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
1074: . N - number of columns (local and global)
1075: . nd - number of block diagonals (global) (optional)
1076: . bs - each element of a diagonal is an bs x bs dense matrix
1077: . diag - optional array of block diagonal numbers (length nd).
1078: For a matrix element A[i,j], where i=row and j=column, the
1079: diagonal number is
1080: $ diag = i/bs - j/bs (integer division)
1081: Set diag=PETSC_NULL on input for PETSc to dynamically allocate memory as
1082: needed (expensive).
1083: - diagv - pointer to actual diagonals (in same order as diag array),
1084: if allocated by user. Otherwise, set diagv=PETSC_NULL on input for PETSc
1085: to control memory allocation.
1087: Output Parameter:
1088: . A - the matrix
1090: Options Database Keys:
1091: . -mat_block_size <bs> - Sets blocksize
1092: . -mat_bdiag_diags <s1,s2,s3,...> - Sets diagonal numbers
1094: Notes:
1095: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
1096: than it must be used on all processors that share the object for that argument.
1098: The parallel matrix is partitioned across the processors by rows, where
1099: each local rectangular matrix is stored in the uniprocessor block
1100: diagonal format. See the users manual for further details.
1102: The user MUST specify either the local or global numbers of rows
1103: (possibly both).
1105: The case bs=1 (conventional diagonal storage) is implemented as
1106: a special case.
1108: Fortran Notes:
1109: Fortran programmers cannot set diagv; this variable is ignored.
1111: Level: intermediate
1113: .keywords: matrix, block, diagonal, parallel, sparse
1115: .seealso: MatCreate(), MatCreateSeqBDiag(), MatSetValues()
1116: @*/
1117: PetscErrorCode MatCreateMPIBDiag(MPI_Comm comm,PetscInt m,PetscInt M,PetscInt N,PetscInt nd,PetscInt bs,const PetscInt diag[],PetscScalar *diagv[],Mat *A)
1118: {
1120: PetscMPIInt size;
1123: MatCreate(comm,A);
1124: MatSetSizes(*A,m,m,M,N);
1125: MPI_Comm_size(comm,&size);
1126: if (size > 1) {
1127: MatSetType(*A,MATMPIBDIAG);
1128: MatMPIBDiagSetPreallocation(*A,nd,bs,diag,diagv);
1129: } else {
1130: MatSetType(*A,MATSEQBDIAG);
1131: MatSeqBDiagSetPreallocation(*A,nd,bs,diag,diagv);
1132: }
1133: return(0);
1134: }
1138: /*@C
1139: MatBDiagGetData - Gets the data for the block diagonal matrix format.
1140: For the parallel case, this returns information for the local submatrix.
1142: Input Parameters:
1143: . mat - the matrix, stored in block diagonal format.
1145: Not Collective
1147: Output Parameters:
1148: + m - number of rows
1149: . n - number of columns
1150: . nd - number of block diagonals
1151: . bs - each element of a diagonal is an bs x bs dense matrix
1152: . bdlen - array of total block lengths of block diagonals
1153: . diag - optional array of block diagonal numbers (length nd).
1154: For a matrix element A[i,j], where i=row and j=column, the
1155: diagonal number is
1156: $ diag = i/bs - j/bs (integer division)
1157: Set diag=PETSC_NULL on input for PETSc to dynamically allocate memory as
1158: needed (expensive).
1159: - diagv - pointer to actual diagonals (in same order as diag array),
1161: Level: advanced
1163: Notes:
1164: See the users manual for further details regarding this storage format.
1166: .keywords: matrix, block, diagonal, get, data
1168: .seealso: MatCreateSeqBDiag(), MatCreateMPIBDiag()
1169: @*/
1170: PetscErrorCode MatBDiagGetData(Mat mat,PetscInt *nd,PetscInt *bs,PetscInt *diag[],PetscInt *bdlen[],PetscScalar ***diagv)
1171: {
1172: Mat_MPIBDiag *pdmat;
1173: Mat_SeqBDiag *dmat = 0;
1174: PetscTruth isseq,ismpi;
1179: PetscTypeCompare((PetscObject)mat,MATSEQBDIAG,&isseq);
1180: PetscTypeCompare((PetscObject)mat,MATMPIBDIAG,&ismpi);
1181: if (isseq) {
1182: dmat = (Mat_SeqBDiag*)mat->data;
1183: } else if (ismpi) {
1184: pdmat = (Mat_MPIBDiag*)mat->data;
1185: dmat = (Mat_SeqBDiag*)pdmat->A->data;
1186: } else SETERRQ(PETSC_ERR_SUP,"Valid only for MATSEQBDIAG and MATMPIBDIAG formats");
1187: *nd = dmat->nd;
1188: *bs = mat->rmap.bs;
1189: *diag = dmat->diag;
1190: *bdlen = dmat->bdlen;
1191: *diagv = dmat->diagv;
1192: return(0);
1193: }
1195: #include petscsys.h
1199: PetscErrorCode MatLoad_MPIBDiag(PetscViewer viewer, MatType type,Mat *newmat)
1200: {
1201: Mat A;
1202: PetscScalar *vals,*svals;
1203: MPI_Comm comm = ((PetscObject)viewer)->comm;
1204: MPI_Status status;
1206: int fd;
1207: PetscMPIInt tag = ((PetscObject)viewer)->tag,rank,size,*sndcounts = 0,*rowners,maxnz,mm;
1208: PetscInt bs,i,nz,j,rstart,rend,*cols;
1209: PetscInt header[4],*rowlengths = 0,M,N,m,Mbs;
1210: PetscInt *ourlens,*procsnz = 0,jj,*mycols,*smycols;
1211: PetscInt extra_rows;
1214: MPI_Comm_size(comm,&size);
1215: MPI_Comm_rank(comm,&rank);
1216: if (!rank) {
1217: PetscViewerBinaryGetDescriptor(viewer,&fd);
1218: PetscBinaryRead(fd,(char *)header,4,PETSC_INT);
1219: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
1220: if (header[3] < 0) {
1221: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format,cannot load as MPIBDiag");
1222: }
1223: }
1224: MPI_Bcast(header+1,3,MPIU_INT,0,comm);
1225: M = header[1]; N = header[2];
1227: bs = 1; /* uses a block size of 1 by default; */
1228: PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);
1230: /*
1231: This code adds extra rows to make sure the number of rows is
1232: divisible by the blocksize
1233: */
1234: Mbs = M/bs;
1235: extra_rows = bs - M + bs*(Mbs);
1236: if (extra_rows == bs) extra_rows = 0;
1237: else Mbs++;
1238: if (extra_rows && !rank) {
1239: PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
1240: }
1242: /* determine ownership of all rows */
1243: m = bs*(Mbs/size + ((Mbs % size) > rank));
1244: PetscMalloc((size+2)*sizeof(PetscInt),&rowners);
1245: mm = (PetscMPIInt)m;
1246: MPI_Allgather(&mm,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1247: rowners[0] = 0;
1248: for (i=2; i<=size; i++) {
1249: rowners[i] += rowners[i-1];
1250: }
1251: rstart = rowners[rank];
1252: rend = rowners[rank+1];
1254: /* distribute row lengths to all processors */
1255: PetscMalloc((rend-rstart)*sizeof(PetscInt),&ourlens);
1256: if (!rank) {
1257: PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);
1258: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
1259: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
1260: PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);
1261: for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1262: MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1263: PetscFree(sndcounts);
1264: } else {
1265: MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1266: }
1268: if (!rank) {
1269: /* calculate the number of nonzeros on each processor */
1270: PetscMalloc(size*sizeof(PetscInt),&procsnz);
1271: PetscMemzero(procsnz,size*sizeof(PetscInt));
1272: for (i=0; i<size; i++) {
1273: for (j=rowners[i]; j<rowners[i+1]; j++) {
1274: procsnz[i] += rowlengths[j];
1275: }
1276: }
1277: PetscFree(rowlengths);
1279: /* determine max buffer needed and allocate it */
1280: maxnz = 0;
1281: for (i=0; i<size; i++) {
1282: maxnz = PetscMax(maxnz,procsnz[i]);
1283: }
1284: PetscMalloc(maxnz*sizeof(PetscInt),&cols);
1286: /* read in my part of the matrix column indices */
1287: nz = procsnz[0];
1288: PetscMalloc(nz*sizeof(PetscInt),&mycols);
1289: if (size == 1) nz -= extra_rows;
1290: PetscBinaryRead(fd,mycols,nz,PETSC_INT);
1291: if (size == 1) for (i=0; i<extra_rows; i++) { mycols[nz+i] = M+i; }
1293: /* read in every one elses and ship off */
1294: for (i=1; i<size-1; i++) {
1295: nz = procsnz[i];
1296: PetscBinaryRead(fd,cols,nz,PETSC_INT);
1297: MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
1298: }
1299: /* read in the stuff for the last proc */
1300: if (size != 1) {
1301: nz = procsnz[size-1] - extra_rows; /* the extra rows are not on the disk */
1302: PetscBinaryRead(fd,cols,nz,PETSC_INT);
1303: for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
1304: MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);
1305: }
1306: PetscFree(cols);
1307: } else {
1308: /* determine buffer space needed for message */
1309: nz = 0;
1310: for (i=0; i<m; i++) {
1311: nz += ourlens[i];
1312: }
1313: PetscMalloc(nz*sizeof(PetscInt),&mycols);
1315: /* receive message of column indices*/
1316: MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
1317: MPI_Get_count(&status,MPIU_INT,&maxnz);
1318: if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
1319: }
1321: MatCreate(comm,newmat);
1322: MatSetSizes(*newmat,m,m,M+extra_rows,N+extra_rows);
1323: MatSetType(*newmat,type);
1324: MatMPIBDiagSetPreallocation(*newmat,0,bs,PETSC_NULL,PETSC_NULL);
1325: A = *newmat;
1327: if (!rank) {
1328: PetscMalloc(maxnz*sizeof(PetscScalar),&vals);
1330: /* read in my part of the matrix numerical values */
1331: nz = procsnz[0];
1332: if (size == 1) nz -= extra_rows;
1333: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1334: if (size == 1) for (i=0; i<extra_rows; i++) { vals[nz+i] = 1.0; }
1336: /* insert into matrix */
1337: jj = rstart;
1338: smycols = mycols;
1339: svals = vals;
1340: for (i=0; i<m; i++) {
1341: MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1342: smycols += ourlens[i];
1343: svals += ourlens[i];
1344: jj++;
1345: }
1347: /* read in other processors (except the last one) and ship out */
1348: for (i=1; i<size-1; i++) {
1349: nz = procsnz[i];
1350: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1351: MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)A)->tag,comm);
1352: }
1353: /* the last proc */
1354: if (size != 1){
1355: nz = procsnz[i] - extra_rows;
1356: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1357: for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
1358: MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,((PetscObject)A)->tag,comm);
1359: }
1360: PetscFree(procsnz);
1361: } else {
1362: /* receive numeric values */
1363: PetscMalloc(nz*sizeof(PetscScalar),&vals);
1365: /* receive message of values*/
1366: MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)A)->tag,comm,&status);
1367: MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
1368: if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
1370: /* insert into matrix */
1371: jj = rstart;
1372: smycols = mycols;
1373: svals = vals;
1374: for (i=0; i<m; i++) {
1375: MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1376: smycols += ourlens[i];
1377: svals += ourlens[i];
1378: jj++;
1379: }
1380: }
1381: PetscFree(ourlens);
1382: PetscFree(vals);
1383: PetscFree(mycols);
1384: PetscFree(rowners);
1386: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1387: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1388: return(0);
1389: }