Actual source code: mpidense.c
1: #define PETSCMAT_DLL
3: /*
4: Basic functions for basic parallel dense matrices.
5: */
6:
7: #include src/mat/impls/dense/mpi/mpidense.h
11: PetscErrorCode MatLUFactorSymbolic_MPIDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
12: {
16: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,fact);
17: return(0);
18: }
22: PetscErrorCode MatCholeskyFactorSymbolic_MPIDense(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
23: {
27: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,fact);
28: return(0);
29: }
33: /*@
35: MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
36: matrix that represents the operator. For sequential matrices it returns itself.
38: Input Parameter:
39: . A - the Seq or MPI dense matrix
41: Output Parameter:
42: . B - the inner matrix
44: Level: intermediate
46: @*/
47: PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
48: {
49: Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
51: PetscTruth flg;
52: PetscMPIInt size;
55: PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);
56: if (!flg) { /* this check sucks! */
57: PetscTypeCompare((PetscObject)A,MATDENSE,&flg);
58: if (flg) {
59: MPI_Comm_size(A->comm,&size);
60: if (size == 1) flg = PETSC_FALSE;
61: }
62: }
63: if (flg) {
64: *B = mat->A;
65: } else {
66: *B = A;
67: }
68: return(0);
69: }
73: PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
74: {
75: Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
77: PetscInt lrow,rstart = A->rmap.rstart,rend = A->rmap.rend;
80: if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
81: lrow = row - rstart;
82: MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);
83: return(0);
84: }
88: PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
89: {
93: if (idx) {PetscFree(*idx);}
94: if (v) {PetscFree(*v);}
95: return(0);
96: }
101: PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
102: {
103: Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
105: PetscInt m = A->rmap.n,rstart = A->rmap.rstart;
106: PetscScalar *array;
107: MPI_Comm comm;
110: if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
112: /* The reuse aspect is not implemented efficiently */
113: if (reuse) { MatDestroy(*B);}
115: PetscObjectGetComm((PetscObject)(mdn->A),&comm);
116: MatGetArray(mdn->A,&array);
117: MatCreate(comm,B);
118: MatSetSizes(*B,m,m,m,m);
119: MatSetType(*B,mdn->A->type_name);
120: MatSeqDenseSetPreallocation(*B,array+m*rstart);
121: MatRestoreArray(mdn->A,&array);
122: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
123: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
124:
125: *iscopy = PETSC_TRUE;
126: return(0);
127: }
132: PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
133: {
134: Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
136: PetscInt i,j,rstart = mat->rmap.rstart,rend = mat->rmap.rend,row;
137: PetscTruth roworiented = A->roworiented;
140: for (i=0; i<m; i++) {
141: if (idxm[i] < 0) continue;
142: if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
143: if (idxm[i] >= rstart && idxm[i] < rend) {
144: row = idxm[i] - rstart;
145: if (roworiented) {
146: MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);
147: } else {
148: for (j=0; j<n; j++) {
149: if (idxn[j] < 0) continue;
150: if (idxn[j] >= mat->cmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
151: MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);
152: }
153: }
154: } else {
155: if (!A->donotstash) {
156: if (roworiented) {
157: MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);
158: } else {
159: MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);
160: }
161: }
162: }
163: }
164: return(0);
165: }
169: PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
170: {
171: Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
173: PetscInt i,j,rstart = mat->rmap.rstart,rend = mat->rmap.rend,row;
176: for (i=0; i<m; i++) {
177: if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
178: if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
179: if (idxm[i] >= rstart && idxm[i] < rend) {
180: row = idxm[i] - rstart;
181: for (j=0; j<n; j++) {
182: if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
183: if (idxn[j] >= mat->cmap.N) {
184: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
185: }
186: MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);
187: }
188: } else {
189: SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
190: }
191: }
192: return(0);
193: }
197: PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
198: {
199: Mat_MPIDense *a = (Mat_MPIDense*)A->data;
203: MatGetArray(a->A,array);
204: return(0);
205: }
209: static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
210: {
211: Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
212: Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
214: PetscInt i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
215: PetscScalar *av,*bv,*v = lmat->v;
216: Mat newmat;
219: ISGetIndices(isrow,&irow);
220: ISGetIndices(iscol,&icol);
221: ISGetLocalSize(isrow,&nrows);
222: ISGetLocalSize(iscol,&ncols);
224: /* No parallel redistribution currently supported! Should really check each index set
225: to comfirm that it is OK. ... Currently supports only submatrix same partitioning as
226: original matrix! */
228: MatGetLocalSize(A,&nlrows,&nlcols);
229: MatGetOwnershipRange(A,&rstart,&rend);
230:
231: /* Check submatrix call */
232: if (scall == MAT_REUSE_MATRIX) {
233: /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
234: /* Really need to test rows and column sizes! */
235: newmat = *B;
236: } else {
237: /* Create and fill new matrix */
238: MatCreate(A->comm,&newmat);
239: MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);
240: MatSetType(newmat,A->type_name);
241: MatMPIDenseSetPreallocation(newmat,PETSC_NULL);
242: }
244: /* Now extract the data pointers and do the copy, column at a time */
245: newmatd = (Mat_MPIDense*)newmat->data;
246: bv = ((Mat_SeqDense *)newmatd->A->data)->v;
247:
248: for (i=0; i<ncols; i++) {
249: av = v + nlrows*icol[i];
250: for (j=0; j<nrows; j++) {
251: *bv++ = av[irow[j] - rstart];
252: }
253: }
255: /* Assemble the matrices so that the correct flags are set */
256: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
257: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
259: /* Free work space */
260: ISRestoreIndices(isrow,&irow);
261: ISRestoreIndices(iscol,&icol);
262: *B = newmat;
263: return(0);
264: }
268: PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
269: {
271: return(0);
272: }
276: PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
277: {
278: Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
279: MPI_Comm comm = mat->comm;
281: PetscInt nstash,reallocs;
282: InsertMode addv;
285: /* make sure all processors are either in INSERTMODE or ADDMODE */
286: MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
287: if (addv == (ADD_VALUES|INSERT_VALUES)) {
288: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
289: }
290: mat->insertmode = addv; /* in case this processor had no cache */
292: MatStashScatterBegin_Private(&mat->stash,mat->rmap.range);
293: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
294: PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
295: return(0);
296: }
300: PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
301: {
302: Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
303: PetscErrorCode ierr;
304: PetscInt i,*row,*col,flg,j,rstart,ncols;
305: PetscMPIInt n;
306: PetscScalar *val;
307: InsertMode addv=mat->insertmode;
310: /* wait on receives */
311: while (1) {
312: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
313: if (!flg) break;
314:
315: for (i=0; i<n;) {
316: /* Now identify the consecutive vals belonging to the same row */
317: for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
318: if (j < n) ncols = j-i;
319: else ncols = n-i;
320: /* Now assemble all these values with a single function call */
321: MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);
322: i = j;
323: }
324: }
325: MatStashScatterEnd_Private(&mat->stash);
326:
327: MatAssemblyBegin(mdn->A,mode);
328: MatAssemblyEnd(mdn->A,mode);
330: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
331: MatSetUpMultiply_MPIDense(mat);
332: }
333: return(0);
334: }
338: PetscErrorCode MatZeroEntries_MPIDense(Mat A)
339: {
341: Mat_MPIDense *l = (Mat_MPIDense*)A->data;
344: MatZeroEntries(l->A);
345: return(0);
346: }
348: /* the code does not do the diagonal entries correctly unless the
349: matrix is square and the column and row owerships are identical.
350: This is a BUG. The only way to fix it seems to be to access
351: mdn->A and mdn->B directly and not through the MatZeroRows()
352: routine.
353: */
356: PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
357: {
358: Mat_MPIDense *l = (Mat_MPIDense*)A->data;
360: PetscInt i,*owners = A->rmap.range;
361: PetscInt *nprocs,j,idx,nsends;
362: PetscInt nmax,*svalues,*starts,*owner,nrecvs;
363: PetscInt *rvalues,tag = A->tag,count,base,slen,*source;
364: PetscInt *lens,*lrows,*values;
365: PetscMPIInt n,imdex,rank = l->rank,size = l->size;
366: MPI_Comm comm = A->comm;
367: MPI_Request *send_waits,*recv_waits;
368: MPI_Status recv_status,*send_status;
369: PetscTruth found;
372: /* first count number of contributors to each processor */
373: PetscMalloc(2*size*sizeof(PetscInt),&nprocs);
374: PetscMemzero(nprocs,2*size*sizeof(PetscInt));
375: PetscMalloc((N+1)*sizeof(PetscInt),&owner); /* see note*/
376: for (i=0; i<N; i++) {
377: idx = rows[i];
378: found = PETSC_FALSE;
379: for (j=0; j<size; j++) {
380: if (idx >= owners[j] && idx < owners[j+1]) {
381: nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
382: }
383: }
384: if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
385: }
386: nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
388: /* inform other processors of number of messages and max length*/
389: PetscMaxSum(comm,nprocs,&nmax,&nrecvs);
391: /* post receives: */
392: PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);
393: PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);
394: for (i=0; i<nrecvs; i++) {
395: MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
396: }
398: /* do sends:
399: 1) starts[i] gives the starting index in svalues for stuff going to
400: the ith processor
401: */
402: PetscMalloc((N+1)*sizeof(PetscInt),&svalues);
403: PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);
404: PetscMalloc((size+1)*sizeof(PetscInt),&starts);
405: starts[0] = 0;
406: for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
407: for (i=0; i<N; i++) {
408: svalues[starts[owner[i]]++] = rows[i];
409: }
411: starts[0] = 0;
412: for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
413: count = 0;
414: for (i=0; i<size; i++) {
415: if (nprocs[2*i+1]) {
416: MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);
417: }
418: }
419: PetscFree(starts);
421: base = owners[rank];
423: /* wait on receives */
424: PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);
425: source = lens + nrecvs;
426: count = nrecvs; slen = 0;
427: while (count) {
428: MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
429: /* unpack receives into our local space */
430: MPI_Get_count(&recv_status,MPIU_INT,&n);
431: source[imdex] = recv_status.MPI_SOURCE;
432: lens[imdex] = n;
433: slen += n;
434: count--;
435: }
436: PetscFree(recv_waits);
437:
438: /* move the data into the send scatter */
439: PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);
440: count = 0;
441: for (i=0; i<nrecvs; i++) {
442: values = rvalues + i*nmax;
443: for (j=0; j<lens[i]; j++) {
444: lrows[count++] = values[j] - base;
445: }
446: }
447: PetscFree(rvalues);
448: PetscFree(lens);
449: PetscFree(owner);
450: PetscFree(nprocs);
451:
452: /* actually zap the local rows */
453: MatZeroRows(l->A,slen,lrows,diag);
454: PetscFree(lrows);
456: /* wait on sends */
457: if (nsends) {
458: PetscMalloc(nsends*sizeof(MPI_Status),&send_status);
459: MPI_Waitall(nsends,send_waits,send_status);
460: PetscFree(send_status);
461: }
462: PetscFree(send_waits);
463: PetscFree(svalues);
465: return(0);
466: }
470: PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
471: {
472: Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
476: VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);
477: VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);
478: MatMult_SeqDense(mdn->A,mdn->lvec,yy);
479: return(0);
480: }
484: PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
485: {
486: Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
490: VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);
491: VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);
492: MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);
493: return(0);
494: }
498: PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
499: {
500: Mat_MPIDense *a = (Mat_MPIDense*)A->data;
502: PetscScalar zero = 0.0;
505: VecSet(yy,zero);
506: MatMultTranspose_SeqDense(a->A,xx,a->lvec);
507: VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);
508: VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);
509: return(0);
510: }
514: PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
515: {
516: Mat_MPIDense *a = (Mat_MPIDense*)A->data;
520: VecCopy(yy,zz);
521: MatMultTranspose_SeqDense(a->A,xx,a->lvec);
522: VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);
523: VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);
524: return(0);
525: }
529: PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
530: {
531: Mat_MPIDense *a = (Mat_MPIDense*)A->data;
532: Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
534: PetscInt len,i,n,m = A->rmap.n,radd;
535: PetscScalar *x,zero = 0.0;
536:
538: VecSet(v,zero);
539: VecGetArray(v,&x);
540: VecGetSize(v,&n);
541: if (n != A->rmap.N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
542: len = PetscMin(a->A->rmap.n,a->A->cmap.n);
543: radd = A->rmap.rstart*m;
544: for (i=0; i<len; i++) {
545: x[i] = aloc->v[radd + i*m + i];
546: }
547: VecRestoreArray(v,&x);
548: return(0);
549: }
553: PetscErrorCode MatDestroy_MPIDense(Mat mat)
554: {
555: Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
560: #if defined(PETSC_USE_LOG)
561: PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap.N,mat->cmap.N);
562: #endif
563: MatStashDestroy_Private(&mat->stash);
564: MatDestroy(mdn->A);
565: if (mdn->lvec) {VecDestroy(mdn->lvec);}
566: if (mdn->Mvctx) {VecScatterDestroy(mdn->Mvctx);}
567: if (mdn->factor) {
568: PetscFree(mdn->factor->temp);
569: PetscFree(mdn->factor->tag);
570: PetscFree(mdn->factor->pivots);
571: PetscFree(mdn->factor);
572: }
573: PetscFree(mdn);
574: PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);
575: PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);
576: return(0);
577: }
581: static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
582: {
583: Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
587: if (mdn->size == 1) {
588: MatView(mdn->A,viewer);
589: }
590: else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
591: return(0);
592: }
596: static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
597: {
598: Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
599: PetscErrorCode ierr;
600: PetscMPIInt size = mdn->size,rank = mdn->rank;
601: PetscViewerType vtype;
602: PetscTruth iascii,isdraw;
603: PetscViewer sviewer;
604: PetscViewerFormat format;
607: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
608: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
609: if (iascii) {
610: PetscViewerGetType(viewer,&vtype);
611: PetscViewerGetFormat(viewer,&format);
612: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
613: MatInfo info;
614: MatGetInfo(mat,MAT_LOCAL,&info);
615: PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap.n,
616: (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);
617: PetscViewerFlush(viewer);
618: VecScatterView(mdn->Mvctx,viewer);
619: return(0);
620: } else if (format == PETSC_VIEWER_ASCII_INFO) {
621: return(0);
622: }
623: } else if (isdraw) {
624: PetscDraw draw;
625: PetscTruth isnull;
627: PetscViewerDrawGetDraw(viewer,0,&draw);
628: PetscDrawIsNull(draw,&isnull);
629: if (isnull) return(0);
630: }
632: if (size == 1) {
633: MatView(mdn->A,viewer);
634: } else {
635: /* assemble the entire matrix onto first processor. */
636: Mat A;
637: PetscInt M = mat->rmap.N,N = mat->cmap.N,m,row,i,nz;
638: PetscInt *cols;
639: PetscScalar *vals;
641: MatCreate(mat->comm,&A);
642: if (!rank) {
643: MatSetSizes(A,M,N,M,N);
644: } else {
645: MatSetSizes(A,0,0,M,N);
646: }
647: /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */
648: MatSetType(A,MATMPIDENSE);
649: MatMPIDenseSetPreallocation(A,PETSC_NULL);
650: PetscLogObjectParent(mat,A);
652: /* Copy the matrix ... This isn't the most efficient means,
653: but it's quick for now */
654: A->insertmode = INSERT_VALUES;
655: row = mat->rmap.rstart; m = mdn->A->rmap.n;
656: for (i=0; i<m; i++) {
657: MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);
658: MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);
659: MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);
660: row++;
661: }
663: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
664: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
665: PetscViewerGetSingleton(viewer,&sviewer);
666: if (!rank) {
667: MatView(((Mat_MPIDense*)(A->data))->A,sviewer);
668: }
669: PetscViewerRestoreSingleton(viewer,&sviewer);
670: PetscViewerFlush(viewer);
671: MatDestroy(A);
672: }
673: return(0);
674: }
678: PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
679: {
681: PetscTruth iascii,isbinary,isdraw,issocket;
682:
684:
685: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
686: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
687: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
688: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
690: if (iascii || issocket || isdraw) {
691: MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);
692: } else if (isbinary) {
693: MatView_MPIDense_Binary(mat,viewer);
694: } else {
695: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
696: }
697: return(0);
698: }
702: PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
703: {
704: Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
705: Mat mdn = mat->A;
707: PetscReal isend[5],irecv[5];
710: info->rows_global = (double)A->rmap.N;
711: info->columns_global = (double)A->cmap.N;
712: info->rows_local = (double)A->rmap.n;
713: info->columns_local = (double)A->cmap.N;
714: info->block_size = 1.0;
715: MatGetInfo(mdn,MAT_LOCAL,info);
716: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
717: isend[3] = info->memory; isend[4] = info->mallocs;
718: if (flag == MAT_LOCAL) {
719: info->nz_used = isend[0];
720: info->nz_allocated = isend[1];
721: info->nz_unneeded = isend[2];
722: info->memory = isend[3];
723: info->mallocs = isend[4];
724: } else if (flag == MAT_GLOBAL_MAX) {
725: MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);
726: info->nz_used = irecv[0];
727: info->nz_allocated = irecv[1];
728: info->nz_unneeded = irecv[2];
729: info->memory = irecv[3];
730: info->mallocs = irecv[4];
731: } else if (flag == MAT_GLOBAL_SUM) {
732: MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);
733: info->nz_used = irecv[0];
734: info->nz_allocated = irecv[1];
735: info->nz_unneeded = irecv[2];
736: info->memory = irecv[3];
737: info->mallocs = irecv[4];
738: }
739: info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */
740: info->fill_ratio_needed = 0;
741: info->factor_mallocs = 0;
742: return(0);
743: }
747: PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op)
748: {
749: Mat_MPIDense *a = (Mat_MPIDense*)A->data;
753: switch (op) {
754: case MAT_NO_NEW_NONZERO_LOCATIONS:
755: case MAT_YES_NEW_NONZERO_LOCATIONS:
756: case MAT_NEW_NONZERO_LOCATION_ERR:
757: case MAT_NEW_NONZERO_ALLOCATION_ERR:
758: case MAT_COLUMNS_SORTED:
759: case MAT_COLUMNS_UNSORTED:
760: MatSetOption(a->A,op);
761: break;
762: case MAT_ROW_ORIENTED:
763: a->roworiented = PETSC_TRUE;
764: MatSetOption(a->A,op);
765: break;
766: case MAT_ROWS_SORTED:
767: case MAT_ROWS_UNSORTED:
768: case MAT_YES_NEW_DIAGONALS:
769: case MAT_USE_HASH_TABLE:
770: PetscInfo(A,"Option ignored\n");
771: break;
772: case MAT_COLUMN_ORIENTED:
773: a->roworiented = PETSC_FALSE;
774: MatSetOption(a->A,op);
775: break;
776: case MAT_IGNORE_OFF_PROC_ENTRIES:
777: a->donotstash = PETSC_TRUE;
778: break;
779: case MAT_NO_NEW_DIAGONALS:
780: SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
781: case MAT_SYMMETRIC:
782: case MAT_STRUCTURALLY_SYMMETRIC:
783: case MAT_NOT_SYMMETRIC:
784: case MAT_NOT_STRUCTURALLY_SYMMETRIC:
785: case MAT_HERMITIAN:
786: case MAT_NOT_HERMITIAN:
787: case MAT_SYMMETRY_ETERNAL:
788: case MAT_NOT_SYMMETRY_ETERNAL:
789: break;
790: default:
791: SETERRQ(PETSC_ERR_SUP,"unknown option");
792: }
793: return(0);
794: }
799: PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
800: {
801: Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
802: Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
803: PetscScalar *l,*r,x,*v;
805: PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->rmap.n,n=mdn->A->cmap.n;
808: MatGetLocalSize(A,&s2,&s3);
809: if (ll) {
810: VecGetLocalSize(ll,&s2a);
811: if (s2a != s2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
812: VecGetArray(ll,&l);
813: for (i=0; i<m; i++) {
814: x = l[i];
815: v = mat->v + i;
816: for (j=0; j<n; j++) { (*v) *= x; v+= m;}
817: }
818: VecRestoreArray(ll,&l);
819: PetscLogFlops(n*m);
820: }
821: if (rr) {
822: VecGetLocalSize(rr,&s3a);
823: if (s3a != s3) SETERRQ2(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
824: VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);
825: VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);
826: VecGetArray(mdn->lvec,&r);
827: for (i=0; i<n; i++) {
828: x = r[i];
829: v = mat->v + i*m;
830: for (j=0; j<m; j++) { (*v++) *= x;}
831: }
832: VecRestoreArray(mdn->lvec,&r);
833: PetscLogFlops(n*m);
834: }
835: return(0);
836: }
840: PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
841: {
842: Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
843: Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
845: PetscInt i,j;
846: PetscReal sum = 0.0;
847: PetscScalar *v = mat->v;
850: if (mdn->size == 1) {
851: MatNorm(mdn->A,type,nrm);
852: } else {
853: if (type == NORM_FROBENIUS) {
854: for (i=0; i<mdn->A->cmap.n*mdn->A->rmap.n; i++) {
855: #if defined(PETSC_USE_COMPLEX)
856: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
857: #else
858: sum += (*v)*(*v); v++;
859: #endif
860: }
861: MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);
862: *nrm = sqrt(*nrm);
863: PetscLogFlops(2*mdn->A->cmap.n*mdn->A->rmap.n);
864: } else if (type == NORM_1) {
865: PetscReal *tmp,*tmp2;
866: PetscMalloc(2*A->cmap.N*sizeof(PetscReal),&tmp);
867: tmp2 = tmp + A->cmap.N;
868: PetscMemzero(tmp,2*A->cmap.N*sizeof(PetscReal));
869: *nrm = 0.0;
870: v = mat->v;
871: for (j=0; j<mdn->A->cmap.n; j++) {
872: for (i=0; i<mdn->A->rmap.n; i++) {
873: tmp[j] += PetscAbsScalar(*v); v++;
874: }
875: }
876: MPI_Allreduce(tmp,tmp2,A->cmap.N,MPIU_REAL,MPI_SUM,A->comm);
877: for (j=0; j<A->cmap.N; j++) {
878: if (tmp2[j] > *nrm) *nrm = tmp2[j];
879: }
880: PetscFree(tmp);
881: PetscLogFlops(A->cmap.n*A->rmap.n);
882: } else if (type == NORM_INFINITY) { /* max row norm */
883: PetscReal ntemp;
884: MatNorm(mdn->A,type,&ntemp);
885: MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);
886: } else {
887: SETERRQ(PETSC_ERR_SUP,"No support for two norm");
888: }
889: }
890: return(0);
891: }
895: PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout)
896: {
897: Mat_MPIDense *a = (Mat_MPIDense*)A->data;
898: Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
899: Mat B;
900: PetscInt M = A->rmap.N,N = A->cmap.N,m,n,*rwork,rstart = A->rmap.rstart;
902: PetscInt j,i;
903: PetscScalar *v;
906: if (!matout && M != N) {
907: SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
908: }
909: MatCreate(A->comm,&B);
910: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);
911: MatSetType(B,A->type_name);
912: MatMPIDenseSetPreallocation(B,PETSC_NULL);
914: m = a->A->rmap.n; n = a->A->cmap.n; v = Aloc->v;
915: PetscMalloc(n*sizeof(PetscInt),&rwork);
916: for (j=0; j<n; j++) {
917: for (i=0; i<m; i++) rwork[i] = rstart + i;
918: MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);
919: v += m;
920: }
921: PetscFree(rwork);
922: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
923: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
924: if (matout) {
925: *matout = B;
926: } else {
927: MatHeaderCopy(A,B);
928: }
929: return(0);
930: }
932: #include petscblaslapack.h
935: PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha)
936: {
937: Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
938: Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
939: PetscScalar oalpha = alpha;
940: PetscBLASInt one = 1,nz = (PetscBLASInt)inA->rmap.n*inA->cmap.N;
944: BLASscal_(&nz,&oalpha,a->v,&one);
945: PetscLogFlops(nz);
946: return(0);
947: }
949: static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
953: PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
954: {
958: MatMPIDenseSetPreallocation(A,0);
959: return(0);
960: }
962: /* -------------------------------------------------------------------*/
963: static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
964: MatGetRow_MPIDense,
965: MatRestoreRow_MPIDense,
966: MatMult_MPIDense,
967: /* 4*/ MatMultAdd_MPIDense,
968: MatMultTranspose_MPIDense,
969: MatMultTransposeAdd_MPIDense,
970: 0,
971: 0,
972: 0,
973: /*10*/ 0,
974: 0,
975: 0,
976: 0,
977: MatTranspose_MPIDense,
978: /*15*/ MatGetInfo_MPIDense,
979: 0,
980: MatGetDiagonal_MPIDense,
981: MatDiagonalScale_MPIDense,
982: MatNorm_MPIDense,
983: /*20*/ MatAssemblyBegin_MPIDense,
984: MatAssemblyEnd_MPIDense,
985: 0,
986: MatSetOption_MPIDense,
987: MatZeroEntries_MPIDense,
988: /*25*/ MatZeroRows_MPIDense,
989: MatLUFactorSymbolic_MPIDense,
990: 0,
991: MatCholeskyFactorSymbolic_MPIDense,
992: 0,
993: /*30*/ MatSetUpPreallocation_MPIDense,
994: 0,
995: 0,
996: MatGetArray_MPIDense,
997: MatRestoreArray_MPIDense,
998: /*35*/ MatDuplicate_MPIDense,
999: 0,
1000: 0,
1001: 0,
1002: 0,
1003: /*40*/ 0,
1004: MatGetSubMatrices_MPIDense,
1005: 0,
1006: MatGetValues_MPIDense,
1007: 0,
1008: /*45*/ 0,
1009: MatScale_MPIDense,
1010: 0,
1011: 0,
1012: 0,
1013: /*50*/ 0,
1014: 0,
1015: 0,
1016: 0,
1017: 0,
1018: /*55*/ 0,
1019: 0,
1020: 0,
1021: 0,
1022: 0,
1023: /*60*/ MatGetSubMatrix_MPIDense,
1024: MatDestroy_MPIDense,
1025: MatView_MPIDense,
1026: 0,
1027: 0,
1028: /*65*/ 0,
1029: 0,
1030: 0,
1031: 0,
1032: 0,
1033: /*70*/ 0,
1034: 0,
1035: 0,
1036: 0,
1037: 0,
1038: /*75*/ 0,
1039: 0,
1040: 0,
1041: 0,
1042: 0,
1043: /*80*/ 0,
1044: 0,
1045: 0,
1046: 0,
1047: /*84*/ MatLoad_MPIDense,
1048: 0,
1049: 0,
1050: 0,
1051: 0,
1052: 0,
1053: /*90*/ 0,
1054: 0,
1055: 0,
1056: 0,
1057: 0,
1058: /*95*/ 0,
1059: 0,
1060: 0,
1061: 0};
1066: PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1067: {
1068: Mat_MPIDense *a;
1072: mat->preallocated = PETSC_TRUE;
1073: /* Note: For now, when data is specified above, this assumes the user correctly
1074: allocates the local dense storage space. We should add error checking. */
1076: a = (Mat_MPIDense*)mat->data;
1077: MatCreate(PETSC_COMM_SELF,&a->A);
1078: MatSetSizes(a->A,mat->rmap.n,mat->cmap.N,mat->rmap.n,mat->cmap.N);
1079: MatSetType(a->A,MATSEQDENSE);
1080: MatSeqDenseSetPreallocation(a->A,data);
1081: PetscLogObjectParent(mat,a->A);
1082: return(0);
1083: }
1086: /*MC
1087: MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
1089: Options Database Keys:
1090: . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
1092: Level: beginner
1094: .seealso: MatCreateMPIDense
1095: M*/
1100: PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat)
1101: {
1102: Mat_MPIDense *a;
1106: PetscNew(Mat_MPIDense,&a);
1107: mat->data = (void*)a;
1108: PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));
1109: mat->factor = 0;
1110: mat->mapping = 0;
1112: a->factor = 0;
1113: mat->insertmode = NOT_SET_VALUES;
1114: MPI_Comm_rank(mat->comm,&a->rank);
1115: MPI_Comm_size(mat->comm,&a->size);
1117: mat->rmap.bs = mat->cmap.bs = 1;
1118: PetscMapInitialize(mat->comm,&mat->rmap);
1119: PetscMapInitialize(mat->comm,&mat->cmap);
1120: a->nvec = mat->cmap.n;
1122: /* build cache for off array entries formed */
1123: a->donotstash = PETSC_FALSE;
1124: MatStashCreate_Private(mat->comm,1,&mat->stash);
1126: /* stuff used for matrix vector multiply */
1127: a->lvec = 0;
1128: a->Mvctx = 0;
1129: a->roworiented = PETSC_TRUE;
1131: PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1132: "MatGetDiagonalBlock_MPIDense",
1133: MatGetDiagonalBlock_MPIDense);
1134: PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1135: "MatMPIDenseSetPreallocation_MPIDense",
1136: MatMPIDenseSetPreallocation_MPIDense);
1137: return(0);
1138: }
1141: /*MC
1142: MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1144: This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1145: and MATMPIDENSE otherwise.
1147: Options Database Keys:
1148: . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1150: Level: beginner
1152: .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1153: M*/
1158: PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A)
1159: {
1161: PetscMPIInt size;
1164: PetscObjectChangeTypeName((PetscObject)A,MATDENSE);
1165: MPI_Comm_size(A->comm,&size);
1166: if (size == 1) {
1167: MatSetType(A,MATSEQDENSE);
1168: } else {
1169: MatSetType(A,MATMPIDENSE);
1170: }
1171: return(0);
1172: }
1177: /*@C
1178: MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1180: Not collective
1182: Input Parameters:
1183: . A - the matrix
1184: - data - optional location of matrix data. Set data=PETSC_NULL for PETSc
1185: to control all matrix memory allocation.
1187: Notes:
1188: The dense format is fully compatible with standard Fortran 77
1189: storage by columns.
1191: The data input variable is intended primarily for Fortran programmers
1192: who wish to allocate their own matrix memory space. Most users should
1193: set data=PETSC_NULL.
1195: Level: intermediate
1197: .keywords: matrix,dense, parallel
1199: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1200: @*/
1201: PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1202: {
1203: PetscErrorCode ierr,(*f)(Mat,PetscScalar *);
1206: PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);
1207: if (f) {
1208: (*f)(mat,data);
1209: }
1210: return(0);
1211: }
1215: /*@C
1216: MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
1218: Collective on MPI_Comm
1220: Input Parameters:
1221: + comm - MPI communicator
1222: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1223: . n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
1224: . M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1225: . N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1226: - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1227: to control all matrix memory allocation.
1229: Output Parameter:
1230: . A - the matrix
1232: Notes:
1233: The dense format is fully compatible with standard Fortran 77
1234: storage by columns.
1236: The data input variable is intended primarily for Fortran programmers
1237: who wish to allocate their own matrix memory space. Most users should
1238: set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
1240: The user MUST specify either the local or global matrix dimensions
1241: (possibly both).
1243: Level: intermediate
1245: .keywords: matrix,dense, parallel
1247: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1248: @*/
1249: PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
1250: {
1252: PetscMPIInt size;
1255: MatCreate(comm,A);
1256: MatSetSizes(*A,m,n,M,N);
1257: MPI_Comm_size(comm,&size);
1258: if (size > 1) {
1259: MatSetType(*A,MATMPIDENSE);
1260: MatMPIDenseSetPreallocation(*A,data);
1261: } else {
1262: MatSetType(*A,MATSEQDENSE);
1263: MatSeqDenseSetPreallocation(*A,data);
1264: }
1265: return(0);
1266: }
1270: static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
1271: {
1272: Mat mat;
1273: Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
1277: *newmat = 0;
1278: MatCreate(A->comm,&mat);
1279: MatSetSizes(mat,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);
1280: MatSetType(mat,A->type_name);
1281: a = (Mat_MPIDense*)mat->data;
1282: PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));
1283: mat->factor = A->factor;
1284: mat->assembled = PETSC_TRUE;
1285: mat->preallocated = PETSC_TRUE;
1287: mat->rmap.rstart = A->rmap.rstart;
1288: mat->rmap.rend = A->rmap.rend;
1289: a->size = oldmat->size;
1290: a->rank = oldmat->rank;
1291: mat->insertmode = NOT_SET_VALUES;
1292: a->nvec = oldmat->nvec;
1293: a->donotstash = oldmat->donotstash;
1294:
1295: PetscMalloc((a->size+1)*sizeof(PetscInt),&mat->rmap.range);
1296: PetscMalloc((a->size+1)*sizeof(PetscInt),&mat->cmap.range);
1297: PetscMemcpy(mat->rmap.range,A->rmap.range,(a->size+1)*sizeof(PetscInt));
1298: PetscMemcpy(mat->cmap.range,A->cmap.range,(a->size+1)*sizeof(PetscInt));
1299: MatStashCreate_Private(A->comm,1,&mat->stash);
1301: MatSetUpMultiply_MPIDense(mat);
1302: MatDuplicate(oldmat->A,cpvalues,&a->A);
1303: PetscLogObjectParent(mat,a->A);
1304: *newmat = mat;
1305: return(0);
1306: }
1308: #include petscsys.h
1312: PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat)
1313: {
1315: PetscMPIInt rank,size;
1316: PetscInt *rowners,i,m,nz,j;
1317: PetscScalar *array,*vals,*vals_ptr;
1318: MPI_Status status;
1321: MPI_Comm_rank(comm,&rank);
1322: MPI_Comm_size(comm,&size);
1324: /* determine ownership of all rows */
1325: m = M/size + ((M % size) > rank);
1326: PetscMalloc((size+2)*sizeof(PetscInt),&rowners);
1327: MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);
1328: rowners[0] = 0;
1329: for (i=2; i<=size; i++) {
1330: rowners[i] += rowners[i-1];
1331: }
1333: MatCreate(comm,newmat);
1334: MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);
1335: MatSetType(*newmat,type);
1336: MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);
1337: MatGetArray(*newmat,&array);
1339: if (!rank) {
1340: PetscMalloc(m*N*sizeof(PetscScalar),&vals);
1342: /* read in my part of the matrix numerical values */
1343: PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);
1344:
1345: /* insert into matrix-by row (this is why cannot directly read into array */
1346: vals_ptr = vals;
1347: for (i=0; i<m; i++) {
1348: for (j=0; j<N; j++) {
1349: array[i + j*m] = *vals_ptr++;
1350: }
1351: }
1353: /* read in other processors and ship out */
1354: for (i=1; i<size; i++) {
1355: nz = (rowners[i+1] - rowners[i])*N;
1356: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1357: MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);
1358: }
1359: } else {
1360: /* receive numeric values */
1361: PetscMalloc(m*N*sizeof(PetscScalar),&vals);
1363: /* receive message of values*/
1364: MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);
1366: /* insert into matrix-by row (this is why cannot directly read into array */
1367: vals_ptr = vals;
1368: for (i=0; i<m; i++) {
1369: for (j=0; j<N; j++) {
1370: array[i + j*m] = *vals_ptr++;
1371: }
1372: }
1373: }
1374: PetscFree(rowners);
1375: PetscFree(vals);
1376: MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
1377: MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
1378: return(0);
1379: }
1383: PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat)
1384: {
1385: Mat A;
1386: PetscScalar *vals,*svals;
1387: MPI_Comm comm = ((PetscObject)viewer)->comm;
1388: MPI_Status status;
1389: PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
1390: PetscInt header[4],*rowlengths = 0,M,N,*cols;
1391: PetscInt *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
1392: PetscInt i,nz,j,rstart,rend;
1393: int fd;
1397: MPI_Comm_size(comm,&size);
1398: MPI_Comm_rank(comm,&rank);
1399: if (!rank) {
1400: PetscViewerBinaryGetDescriptor(viewer,&fd);
1401: PetscBinaryRead(fd,(char *)header,4,PETSC_INT);
1402: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
1403: }
1405: MPI_Bcast(header+1,3,MPIU_INT,0,comm);
1406: M = header[1]; N = header[2]; nz = header[3];
1408: /*
1409: Handle case where matrix is stored on disk as a dense matrix
1410: */
1411: if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1412: MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);
1413: return(0);
1414: }
1416: /* determine ownership of all rows */
1417: m = M/size + ((M % size) > rank);
1418: PetscMalloc((size+2)*sizeof(PetscInt),&rowners);
1419: MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1420: rowners[0] = 0;
1421: for (i=2; i<=size; i++) {
1422: rowners[i] += rowners[i-1];
1423: }
1424: rstart = rowners[rank];
1425: rend = rowners[rank+1];
1427: /* distribute row lengths to all processors */
1428: PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);
1429: offlens = ourlens + (rend-rstart);
1430: if (!rank) {
1431: PetscMalloc(M*sizeof(PetscInt),&rowlengths);
1432: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
1433: PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);
1434: for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1435: MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1436: PetscFree(sndcounts);
1437: } else {
1438: MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1439: }
1441: if (!rank) {
1442: /* calculate the number of nonzeros on each processor */
1443: PetscMalloc(size*sizeof(PetscInt),&procsnz);
1444: PetscMemzero(procsnz,size*sizeof(PetscInt));
1445: for (i=0; i<size; i++) {
1446: for (j=rowners[i]; j< rowners[i+1]; j++) {
1447: procsnz[i] += rowlengths[j];
1448: }
1449: }
1450: PetscFree(rowlengths);
1452: /* determine max buffer needed and allocate it */
1453: maxnz = 0;
1454: for (i=0; i<size; i++) {
1455: maxnz = PetscMax(maxnz,procsnz[i]);
1456: }
1457: PetscMalloc(maxnz*sizeof(PetscInt),&cols);
1459: /* read in my part of the matrix column indices */
1460: nz = procsnz[0];
1461: PetscMalloc(nz*sizeof(PetscInt),&mycols);
1462: PetscBinaryRead(fd,mycols,nz,PETSC_INT);
1464: /* read in every one elses and ship off */
1465: for (i=1; i<size; i++) {
1466: nz = procsnz[i];
1467: PetscBinaryRead(fd,cols,nz,PETSC_INT);
1468: MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
1469: }
1470: PetscFree(cols);
1471: } else {
1472: /* determine buffer space needed for message */
1473: nz = 0;
1474: for (i=0; i<m; i++) {
1475: nz += ourlens[i];
1476: }
1477: PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);
1479: /* receive message of column indices*/
1480: MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
1481: MPI_Get_count(&status,MPIU_INT,&maxnz);
1482: if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
1483: }
1485: /* loop over local rows, determining number of off diagonal entries */
1486: PetscMemzero(offlens,m*sizeof(PetscInt));
1487: jj = 0;
1488: for (i=0; i<m; i++) {
1489: for (j=0; j<ourlens[i]; j++) {
1490: if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
1491: jj++;
1492: }
1493: }
1495: /* create our matrix */
1496: for (i=0; i<m; i++) {
1497: ourlens[i] -= offlens[i];
1498: }
1499: MatCreate(comm,newmat);
1500: MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);
1501: MatSetType(*newmat,type);
1502: MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);
1503: A = *newmat;
1504: for (i=0; i<m; i++) {
1505: ourlens[i] += offlens[i];
1506: }
1508: if (!rank) {
1509: PetscMalloc(maxnz*sizeof(PetscScalar),&vals);
1511: /* read in my part of the matrix numerical values */
1512: nz = procsnz[0];
1513: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1514:
1515: /* insert into matrix */
1516: jj = rstart;
1517: smycols = mycols;
1518: svals = vals;
1519: for (i=0; i<m; i++) {
1520: MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1521: smycols += ourlens[i];
1522: svals += ourlens[i];
1523: jj++;
1524: }
1526: /* read in other processors and ship out */
1527: for (i=1; i<size; i++) {
1528: nz = procsnz[i];
1529: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1530: MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);
1531: }
1532: PetscFree(procsnz);
1533: } else {
1534: /* receive numeric values */
1535: PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);
1537: /* receive message of values*/
1538: MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);
1539: MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
1540: if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
1542: /* insert into matrix */
1543: jj = rstart;
1544: smycols = mycols;
1545: svals = vals;
1546: for (i=0; i<m; i++) {
1547: MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1548: smycols += ourlens[i];
1549: svals += ourlens[i];
1550: jj++;
1551: }
1552: }
1553: PetscFree(ourlens);
1554: PetscFree(vals);
1555: PetscFree(mycols);
1556: PetscFree(rowners);
1558: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1559: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1560: return(0);
1561: }