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  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:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","",PETSC_NULL);
577:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","",PETSC_NULL);
578:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","",PETSC_NULL);
579:   return(0);
580: }

584: static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
585: {
586:   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;

590:   if (mdn->size == 1) {
591:     MatView(mdn->A,viewer);
592:   }
593:   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
594:   return(0);
595: }

599: static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
600: {
601:   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
602:   PetscErrorCode    ierr;
603:   PetscMPIInt       size = mdn->size,rank = mdn->rank;
604:   PetscViewerType   vtype;
605:   PetscTruth        iascii,isdraw;
606:   PetscViewer       sviewer;
607:   PetscViewerFormat format;

610:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
611:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
612:   if (iascii) {
613:     PetscViewerGetType(viewer,&vtype);
614:     PetscViewerGetFormat(viewer,&format);
615:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
616:       MatInfo info;
617:       MatGetInfo(mat,MAT_LOCAL,&info);
618:       PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap.n,
619:                    (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);
620:       PetscViewerFlush(viewer);
621:       VecScatterView(mdn->Mvctx,viewer);
622:       return(0);
623:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
624:       return(0);
625:     }
626:   } else if (isdraw) {
627:     PetscDraw  draw;
628:     PetscTruth isnull;

630:     PetscViewerDrawGetDraw(viewer,0,&draw);
631:     PetscDrawIsNull(draw,&isnull);
632:     if (isnull) return(0);
633:   }

635:   if (size == 1) {
636:     MatView(mdn->A,viewer);
637:   } else {
638:     /* assemble the entire matrix onto first processor. */
639:     Mat         A;
640:     PetscInt    M = mat->rmap.N,N = mat->cmap.N,m,row,i,nz;
641:     PetscInt    *cols;
642:     PetscScalar *vals;

644:     MatCreate(mat->comm,&A);
645:     if (!rank) {
646:       MatSetSizes(A,M,N,M,N);
647:     } else {
648:       MatSetSizes(A,0,0,M,N);
649:     }
650:     /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */
651:     MatSetType(A,MATMPIDENSE);
652:     MatMPIDenseSetPreallocation(A,PETSC_NULL);
653:     PetscLogObjectParent(mat,A);

655:     /* Copy the matrix ... This isn't the most efficient means,
656:        but it's quick for now */
657:     A->insertmode = INSERT_VALUES;
658:     row = mat->rmap.rstart; m = mdn->A->rmap.n;
659:     for (i=0; i<m; i++) {
660:       MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);
661:       MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);
662:       MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);
663:       row++;
664:     }

666:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
667:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
668:     PetscViewerGetSingleton(viewer,&sviewer);
669:     if (!rank) {
670:       MatView(((Mat_MPIDense*)(A->data))->A,sviewer);
671:     }
672:     PetscViewerRestoreSingleton(viewer,&sviewer);
673:     PetscViewerFlush(viewer);
674:     MatDestroy(A);
675:   }
676:   return(0);
677: }

681: PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
682: {
684:   PetscTruth     iascii,isbinary,isdraw,issocket;
685: 
687: 
688:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
689:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
690:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
691:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);

693:   if (iascii || issocket || isdraw) {
694:     MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);
695:   } else if (isbinary) {
696:     MatView_MPIDense_Binary(mat,viewer);
697:   } else {
698:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
699:   }
700:   return(0);
701: }

705: PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
706: {
707:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
708:   Mat            mdn = mat->A;
710:   PetscReal      isend[5],irecv[5];

713:   info->rows_global    = (double)A->rmap.N;
714:   info->columns_global = (double)A->cmap.N;
715:   info->rows_local     = (double)A->rmap.n;
716:   info->columns_local  = (double)A->cmap.N;
717:   info->block_size     = 1.0;
718:   MatGetInfo(mdn,MAT_LOCAL,info);
719:   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
720:   isend[3] = info->memory;  isend[4] = info->mallocs;
721:   if (flag == MAT_LOCAL) {
722:     info->nz_used      = isend[0];
723:     info->nz_allocated = isend[1];
724:     info->nz_unneeded  = isend[2];
725:     info->memory       = isend[3];
726:     info->mallocs      = isend[4];
727:   } else if (flag == MAT_GLOBAL_MAX) {
728:     MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);
729:     info->nz_used      = irecv[0];
730:     info->nz_allocated = irecv[1];
731:     info->nz_unneeded  = irecv[2];
732:     info->memory       = irecv[3];
733:     info->mallocs      = irecv[4];
734:   } else if (flag == MAT_GLOBAL_SUM) {
735:     MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);
736:     info->nz_used      = irecv[0];
737:     info->nz_allocated = irecv[1];
738:     info->nz_unneeded  = irecv[2];
739:     info->memory       = irecv[3];
740:     info->mallocs      = irecv[4];
741:   }
742:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
743:   info->fill_ratio_needed = 0;
744:   info->factor_mallocs    = 0;
745:   return(0);
746: }

750: PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op)
751: {
752:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

756:   switch (op) {
757:   case MAT_NO_NEW_NONZERO_LOCATIONS:
758:   case MAT_YES_NEW_NONZERO_LOCATIONS:
759:   case MAT_NEW_NONZERO_LOCATION_ERR:
760:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
761:   case MAT_COLUMNS_SORTED:
762:   case MAT_COLUMNS_UNSORTED:
763:     MatSetOption(a->A,op);
764:     break;
765:   case MAT_ROW_ORIENTED:
766:     a->roworiented = PETSC_TRUE;
767:     MatSetOption(a->A,op);
768:     break;
769:   case MAT_ROWS_SORTED:
770:   case MAT_ROWS_UNSORTED:
771:   case MAT_YES_NEW_DIAGONALS:
772:   case MAT_USE_HASH_TABLE:
773:     PetscInfo1(A,"Option %d ignored\n",op);
774:     break;
775:   case MAT_COLUMN_ORIENTED:
776:     a->roworiented = PETSC_FALSE;
777:     MatSetOption(a->A,op);
778:     break;
779:   case MAT_IGNORE_OFF_PROC_ENTRIES:
780:     a->donotstash = PETSC_TRUE;
781:     break;
782:   case MAT_NO_NEW_DIAGONALS:
783:     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
784:   case MAT_SYMMETRIC:
785:   case MAT_STRUCTURALLY_SYMMETRIC:
786:   case MAT_NOT_SYMMETRIC:
787:   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
788:   case MAT_HERMITIAN:
789:   case MAT_NOT_HERMITIAN:
790:   case MAT_SYMMETRY_ETERNAL:
791:   case MAT_NOT_SYMMETRY_ETERNAL:
792:     break;
793:   default:
794:     SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
795:   }
796:   return(0);
797: }


802: PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
803: {
804:   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
805:   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
806:   PetscScalar    *l,*r,x,*v;
808:   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap.n,n=mdn->A->cmap.n;

811:   MatGetLocalSize(A,&s2,&s3);
812:   if (ll) {
813:     VecGetLocalSize(ll,&s2a);
814:     if (s2a != s2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
815:     VecGetArray(ll,&l);
816:     for (i=0; i<m; i++) {
817:       x = l[i];
818:       v = mat->v + i;
819:       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
820:     }
821:     VecRestoreArray(ll,&l);
822:     PetscLogFlops(n*m);
823:   }
824:   if (rr) {
825:     VecGetLocalSize(rr,&s3a);
826:     if (s3a != s3) SETERRQ2(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
827:     VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);
828:     VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);
829:     VecGetArray(mdn->lvec,&r);
830:     for (i=0; i<n; i++) {
831:       x = r[i];
832:       v = mat->v + i*m;
833:       for (j=0; j<m; j++) { (*v++) *= x;}
834:     }
835:     VecRestoreArray(mdn->lvec,&r);
836:     PetscLogFlops(n*m);
837:   }
838:   return(0);
839: }

843: PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
844: {
845:   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
846:   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
848:   PetscInt       i,j;
849:   PetscReal      sum = 0.0;
850:   PetscScalar    *v = mat->v;

853:   if (mdn->size == 1) {
854:      MatNorm(mdn->A,type,nrm);
855:   } else {
856:     if (type == NORM_FROBENIUS) {
857:       for (i=0; i<mdn->A->cmap.n*mdn->A->rmap.n; i++) {
858: #if defined(PETSC_USE_COMPLEX)
859:         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
860: #else
861:         sum += (*v)*(*v); v++;
862: #endif
863:       }
864:       MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);
865:       *nrm = sqrt(*nrm);
866:       PetscLogFlops(2*mdn->A->cmap.n*mdn->A->rmap.n);
867:     } else if (type == NORM_1) {
868:       PetscReal *tmp,*tmp2;
869:       PetscMalloc(2*A->cmap.N*sizeof(PetscReal),&tmp);
870:       tmp2 = tmp + A->cmap.N;
871:       PetscMemzero(tmp,2*A->cmap.N*sizeof(PetscReal));
872:       *nrm = 0.0;
873:       v = mat->v;
874:       for (j=0; j<mdn->A->cmap.n; j++) {
875:         for (i=0; i<mdn->A->rmap.n; i++) {
876:           tmp[j] += PetscAbsScalar(*v);  v++;
877:         }
878:       }
879:       MPI_Allreduce(tmp,tmp2,A->cmap.N,MPIU_REAL,MPI_SUM,A->comm);
880:       for (j=0; j<A->cmap.N; j++) {
881:         if (tmp2[j] > *nrm) *nrm = tmp2[j];
882:       }
883:       PetscFree(tmp);
884:       PetscLogFlops(A->cmap.n*A->rmap.n);
885:     } else if (type == NORM_INFINITY) { /* max row norm */
886:       PetscReal ntemp;
887:       MatNorm(mdn->A,type,&ntemp);
888:       MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);
889:     } else {
890:       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
891:     }
892:   }
893:   return(0);
894: }

898: PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout)
899: {
900:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
901:   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
902:   Mat            B;
903:   PetscInt       M = A->rmap.N,N = A->cmap.N,m,n,*rwork,rstart = A->rmap.rstart;
905:   PetscInt       j,i;
906:   PetscScalar    *v;

909:   if (!matout && M != N) {
910:     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
911:   }
912:   MatCreate(A->comm,&B);
913:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);
914:   MatSetType(B,A->type_name);
915:   MatMPIDenseSetPreallocation(B,PETSC_NULL);

917:   m = a->A->rmap.n; n = a->A->cmap.n; v = Aloc->v;
918:   PetscMalloc(n*sizeof(PetscInt),&rwork);
919:   for (j=0; j<n; j++) {
920:     for (i=0; i<m; i++) rwork[i] = rstart + i;
921:     MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);
922:     v   += m;
923:   }
924:   PetscFree(rwork);
925:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
926:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
927:   if (matout) {
928:     *matout = B;
929:   } else {
930:     MatHeaderCopy(A,B);
931:   }
932:   return(0);
933: }

935:  #include petscblaslapack.h
938: PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha)
939: {
940:   Mat_MPIDense   *A = (Mat_MPIDense*)inA->data;
941:   Mat_SeqDense   *a = (Mat_SeqDense*)A->A->data;
942:   PetscScalar    oalpha = alpha;
943:   PetscBLASInt   one = 1,nz = (PetscBLASInt)inA->rmap.n*inA->cmap.N;

947:   BLASscal_(&nz,&oalpha,a->v,&one);
948:   PetscLogFlops(nz);
949:   return(0);
950: }

952: static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);

956: PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
957: {

961:    MatMPIDenseSetPreallocation(A,0);
962:   return(0);
963: }

967: PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
968: {
970:   PetscInt       m=A->rmap.n,n=B->cmap.n;
971:   Mat            Cmat;

974:   if (A->cmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap.n %d != B->rmap.n %d\n",A->cmap.n,B->rmap.n);
975:   MatCreate(B->comm,&Cmat);
976:   MatSetSizes(Cmat,m,n,A->rmap.N,B->cmap.N);
977:   MatSetType(Cmat,MATMPIDENSE);
978:   MatMPIDenseSetPreallocation(Cmat,PETSC_NULL);
979:   MatAssemblyBegin(Cmat,MAT_FINAL_ASSEMBLY);
980:   MatAssemblyEnd(Cmat,MAT_FINAL_ASSEMBLY);
981:   *C = Cmat;
982:   return(0);
983: }

985: /* -------------------------------------------------------------------*/
986: static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
987:        MatGetRow_MPIDense,
988:        MatRestoreRow_MPIDense,
989:        MatMult_MPIDense,
990: /* 4*/ MatMultAdd_MPIDense,
991:        MatMultTranspose_MPIDense,
992:        MatMultTransposeAdd_MPIDense,
993:        0,
994:        0,
995:        0,
996: /*10*/ 0,
997:        0,
998:        0,
999:        0,
1000:        MatTranspose_MPIDense,
1001: /*15*/ MatGetInfo_MPIDense,
1002:        MatEqual_MPIDense,
1003:        MatGetDiagonal_MPIDense,
1004:        MatDiagonalScale_MPIDense,
1005:        MatNorm_MPIDense,
1006: /*20*/ MatAssemblyBegin_MPIDense,
1007:        MatAssemblyEnd_MPIDense,
1008:        0,
1009:        MatSetOption_MPIDense,
1010:        MatZeroEntries_MPIDense,
1011: /*25*/ MatZeroRows_MPIDense,
1012:        MatLUFactorSymbolic_MPIDense,
1013:        0,
1014:        MatCholeskyFactorSymbolic_MPIDense,
1015:        0,
1016: /*30*/ MatSetUpPreallocation_MPIDense,
1017:        0,
1018:        0,
1019:        MatGetArray_MPIDense,
1020:        MatRestoreArray_MPIDense,
1021: /*35*/ MatDuplicate_MPIDense,
1022:        0,
1023:        0,
1024:        0,
1025:        0,
1026: /*40*/ 0,
1027:        MatGetSubMatrices_MPIDense,
1028:        0,
1029:        MatGetValues_MPIDense,
1030:        0,
1031: /*45*/ 0,
1032:        MatScale_MPIDense,
1033:        0,
1034:        0,
1035:        0,
1036: /*50*/ 0,
1037:        0,
1038:        0,
1039:        0,
1040:        0,
1041: /*55*/ 0,
1042:        0,
1043:        0,
1044:        0,
1045:        0,
1046: /*60*/ MatGetSubMatrix_MPIDense,
1047:        MatDestroy_MPIDense,
1048:        MatView_MPIDense,
1049:        0,
1050:        0,
1051: /*65*/ 0,
1052:        0,
1053:        0,
1054:        0,
1055:        0,
1056: /*70*/ 0,
1057:        0,
1058:        0,
1059:        0,
1060:        0,
1061: /*75*/ 0,
1062:        0,
1063:        0,
1064:        0,
1065:        0,
1066: /*80*/ 0,
1067:        0,
1068:        0,
1069:        0,
1070: /*84*/ MatLoad_MPIDense,
1071:        0,
1072:        0,
1073:        0,
1074:        0,
1075:        0,
1076: /*90*/ 0,
1077:        MatMatMultSymbolic_MPIDense_MPIDense,
1078:        0,
1079:        0,
1080:        0,
1081: /*95*/ 0,
1082:        0,
1083:        0,
1084:        0};

1089: PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1090: {
1091:   Mat_MPIDense   *a;

1095:   mat->preallocated = PETSC_TRUE;
1096:   /* Note:  For now, when data is specified above, this assumes the user correctly
1097:    allocates the local dense storage space.  We should add error checking. */

1099:   a    = (Mat_MPIDense*)mat->data;
1100:   MatCreate(PETSC_COMM_SELF,&a->A);
1101:   MatSetSizes(a->A,mat->rmap.n,mat->cmap.N,mat->rmap.n,mat->cmap.N);
1102:   MatSetType(a->A,MATSEQDENSE);
1103:   MatSeqDenseSetPreallocation(a->A,data);
1104:   PetscLogObjectParent(mat,a->A);
1105:   return(0);
1106: }

1109: /*MC
1110:    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.

1112:    Options Database Keys:
1113: . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()

1115:   Level: beginner

1117: .seealso: MatCreateMPIDense
1118: M*/

1123: PetscErrorCode  MatCreate_MPIDense(Mat mat)
1124: {
1125:   Mat_MPIDense   *a;

1129:   PetscNew(Mat_MPIDense,&a);
1130:   mat->data         = (void*)a;
1131:   PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));
1132:   mat->factor       = 0;
1133:   mat->mapping      = 0;

1135:   a->factor       = 0;
1136:   mat->insertmode = NOT_SET_VALUES;
1137:   MPI_Comm_rank(mat->comm,&a->rank);
1138:   MPI_Comm_size(mat->comm,&a->size);

1140:   mat->rmap.bs = mat->cmap.bs = 1;
1141:   PetscMapInitialize(mat->comm,&mat->rmap);
1142:   PetscMapInitialize(mat->comm,&mat->cmap);
1143:   a->nvec = mat->cmap.n;

1145:   /* build cache for off array entries formed */
1146:   a->donotstash = PETSC_FALSE;
1147:   MatStashCreate_Private(mat->comm,1,&mat->stash);

1149:   /* stuff used for matrix vector multiply */
1150:   a->lvec        = 0;
1151:   a->Mvctx       = 0;
1152:   a->roworiented = PETSC_TRUE;

1154:   PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1155:                                      "MatGetDiagonalBlock_MPIDense",
1156:                                      MatGetDiagonalBlock_MPIDense);
1157:   PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1158:                                      "MatMPIDenseSetPreallocation_MPIDense",
1159:                                      MatMPIDenseSetPreallocation_MPIDense);
1160:   PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",
1161:                                      "MatMatMult_MPIAIJ_MPIDense",
1162:                                       MatMatMult_MPIAIJ_MPIDense);
1163:   PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",
1164:                                      "MatMatMultSymbolic_MPIAIJ_MPIDense",
1165:                                       MatMatMultSymbolic_MPIAIJ_MPIDense);
1166:   PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",
1167:                                      "MatMatMultNumeric_MPIAIJ_MPIDense",
1168:                                       MatMatMultNumeric_MPIAIJ_MPIDense);
1169:   PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);
1170:   return(0);
1171: }

1174: /*MC
1175:    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.

1177:    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1178:    and MATMPIDENSE otherwise.

1180:    Options Database Keys:
1181: . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()

1183:   Level: beginner

1185: .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1186: M*/

1191: PetscErrorCode  MatCreate_Dense(Mat A)
1192: {
1194:   PetscMPIInt    size;

1197:   MPI_Comm_size(A->comm,&size);
1198:   if (size == 1) {
1199:     MatSetType(A,MATSEQDENSE);
1200:   } else {
1201:     MatSetType(A,MATMPIDENSE);
1202:   }
1203:   return(0);
1204: }

1209: /*@C
1210:    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries

1212:    Not collective

1214:    Input Parameters:
1215: .  A - the matrix
1216: -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1217:    to control all matrix memory allocation.

1219:    Notes:
1220:    The dense format is fully compatible with standard Fortran 77
1221:    storage by columns.

1223:    The data input variable is intended primarily for Fortran programmers
1224:    who wish to allocate their own matrix memory space.  Most users should
1225:    set data=PETSC_NULL.

1227:    Level: intermediate

1229: .keywords: matrix,dense, parallel

1231: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1232: @*/
1233: PetscErrorCode  MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1234: {
1235:   PetscErrorCode ierr,(*f)(Mat,PetscScalar *);

1238:   PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);
1239:   if (f) {
1240:     (*f)(mat,data);
1241:   }
1242:   return(0);
1243: }

1247: /*@C
1248:    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.

1250:    Collective on MPI_Comm

1252:    Input Parameters:
1253: +  comm - MPI communicator
1254: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1255: .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
1256: .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1257: .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1258: -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1259:    to control all matrix memory allocation.

1261:    Output Parameter:
1262: .  A - the matrix

1264:    Notes:
1265:    The dense format is fully compatible with standard Fortran 77
1266:    storage by columns.

1268:    The data input variable is intended primarily for Fortran programmers
1269:    who wish to allocate their own matrix memory space.  Most users should
1270:    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).

1272:    The user MUST specify either the local or global matrix dimensions
1273:    (possibly both).

1275:    Level: intermediate

1277: .keywords: matrix,dense, parallel

1279: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1280: @*/
1281: PetscErrorCode  MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
1282: {
1284:   PetscMPIInt    size;

1287:   MatCreate(comm,A);
1288:   MatSetSizes(*A,m,n,M,N);
1289:   MPI_Comm_size(comm,&size);
1290:   if (size > 1) {
1291:     MatSetType(*A,MATMPIDENSE);
1292:     MatMPIDenseSetPreallocation(*A,data);
1293:   } else {
1294:     MatSetType(*A,MATSEQDENSE);
1295:     MatSeqDenseSetPreallocation(*A,data);
1296:   }
1297:   return(0);
1298: }

1302: static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
1303: {
1304:   Mat            mat;
1305:   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;

1309:   *newmat       = 0;
1310:   MatCreate(A->comm,&mat);
1311:   MatSetSizes(mat,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);
1312:   MatSetType(mat,A->type_name);
1313:   a                 = (Mat_MPIDense*)mat->data;
1314:   PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));
1315:   mat->factor       = A->factor;
1316:   mat->assembled    = PETSC_TRUE;
1317:   mat->preallocated = PETSC_TRUE;

1319:   mat->rmap.rstart     = A->rmap.rstart;
1320:   mat->rmap.rend       = A->rmap.rend;
1321:   a->size              = oldmat->size;
1322:   a->rank              = oldmat->rank;
1323:   mat->insertmode      = NOT_SET_VALUES;
1324:   a->nvec              = oldmat->nvec;
1325:   a->donotstash        = oldmat->donotstash;
1326: 
1327:   PetscMemcpy(mat->rmap.range,A->rmap.range,(a->size+1)*sizeof(PetscInt));
1328:   PetscMemcpy(mat->cmap.range,A->cmap.range,(a->size+1)*sizeof(PetscInt));
1329:   MatStashCreate_Private(A->comm,1,&mat->stash);

1331:   MatSetUpMultiply_MPIDense(mat);
1332:   MatDuplicate(oldmat->A,cpvalues,&a->A);
1333:   PetscLogObjectParent(mat,a->A);
1334:   *newmat = mat;
1335:   return(0);
1336: }

1338:  #include petscsys.h

1342: PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat)
1343: {
1345:   PetscMPIInt    rank,size;
1346:   PetscInt       *rowners,i,m,nz,j;
1347:   PetscScalar    *array,*vals,*vals_ptr;
1348:   MPI_Status     status;

1351:   MPI_Comm_rank(comm,&rank);
1352:   MPI_Comm_size(comm,&size);

1354:   /* determine ownership of all rows */
1355:   m          = M/size + ((M % size) > rank);
1356:   PetscMalloc((size+2)*sizeof(PetscInt),&rowners);
1357:   MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);
1358:   rowners[0] = 0;
1359:   for (i=2; i<=size; i++) {
1360:     rowners[i] += rowners[i-1];
1361:   }

1363:   MatCreate(comm,newmat);
1364:   MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);
1365:   MatSetType(*newmat,type);
1366:   MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);
1367:   MatGetArray(*newmat,&array);

1369:   if (!rank) {
1370:     PetscMalloc(m*N*sizeof(PetscScalar),&vals);

1372:     /* read in my part of the matrix numerical values  */
1373:     PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);
1374: 
1375:     /* insert into matrix-by row (this is why cannot directly read into array */
1376:     vals_ptr = vals;
1377:     for (i=0; i<m; i++) {
1378:       for (j=0; j<N; j++) {
1379:         array[i + j*m] = *vals_ptr++;
1380:       }
1381:     }

1383:     /* read in other processors and ship out */
1384:     for (i=1; i<size; i++) {
1385:       nz   = (rowners[i+1] - rowners[i])*N;
1386:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1387:       MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);
1388:     }
1389:   } else {
1390:     /* receive numeric values */
1391:     PetscMalloc(m*N*sizeof(PetscScalar),&vals);

1393:     /* receive message of values*/
1394:     MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);

1396:     /* insert into matrix-by row (this is why cannot directly read into array */
1397:     vals_ptr = vals;
1398:     for (i=0; i<m; i++) {
1399:       for (j=0; j<N; j++) {
1400:         array[i + j*m] = *vals_ptr++;
1401:       }
1402:     }
1403:   }
1404:   PetscFree(rowners);
1405:   PetscFree(vals);
1406:   MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
1407:   MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
1408:   return(0);
1409: }

1413: PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat)
1414: {
1415:   Mat            A;
1416:   PetscScalar    *vals,*svals;
1417:   MPI_Comm       comm = ((PetscObject)viewer)->comm;
1418:   MPI_Status     status;
1419:   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
1420:   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
1421:   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
1422:   PetscInt       i,nz,j,rstart,rend;
1423:   int            fd;

1427:   MPI_Comm_size(comm,&size);
1428:   MPI_Comm_rank(comm,&rank);
1429:   if (!rank) {
1430:     PetscViewerBinaryGetDescriptor(viewer,&fd);
1431:     PetscBinaryRead(fd,(char *)header,4,PETSC_INT);
1432:     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
1433:   }

1435:   MPI_Bcast(header+1,3,MPIU_INT,0,comm);
1436:   M = header[1]; N = header[2]; nz = header[3];

1438:   /*
1439:        Handle case where matrix is stored on disk as a dense matrix 
1440:   */
1441:   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1442:     MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);
1443:     return(0);
1444:   }

1446:   /* determine ownership of all rows */
1447:   m          = M/size + ((M % size) > rank);
1448:   PetscMalloc((size+2)*sizeof(PetscInt),&rowners);
1449:   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1450:   rowners[0] = 0;
1451:   for (i=2; i<=size; i++) {
1452:     rowners[i] += rowners[i-1];
1453:   }
1454:   rstart = rowners[rank];
1455:   rend   = rowners[rank+1];

1457:   /* distribute row lengths to all processors */
1458:   PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);
1459:   offlens = ourlens + (rend-rstart);
1460:   if (!rank) {
1461:     PetscMalloc(M*sizeof(PetscInt),&rowlengths);
1462:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
1463:     PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);
1464:     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1465:     MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1466:     PetscFree(sndcounts);
1467:   } else {
1468:     MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1469:   }

1471:   if (!rank) {
1472:     /* calculate the number of nonzeros on each processor */
1473:     PetscMalloc(size*sizeof(PetscInt),&procsnz);
1474:     PetscMemzero(procsnz,size*sizeof(PetscInt));
1475:     for (i=0; i<size; i++) {
1476:       for (j=rowners[i]; j< rowners[i+1]; j++) {
1477:         procsnz[i] += rowlengths[j];
1478:       }
1479:     }
1480:     PetscFree(rowlengths);

1482:     /* determine max buffer needed and allocate it */
1483:     maxnz = 0;
1484:     for (i=0; i<size; i++) {
1485:       maxnz = PetscMax(maxnz,procsnz[i]);
1486:     }
1487:     PetscMalloc(maxnz*sizeof(PetscInt),&cols);

1489:     /* read in my part of the matrix column indices  */
1490:     nz = procsnz[0];
1491:     PetscMalloc(nz*sizeof(PetscInt),&mycols);
1492:     PetscBinaryRead(fd,mycols,nz,PETSC_INT);

1494:     /* read in every one elses and ship off */
1495:     for (i=1; i<size; i++) {
1496:       nz   = procsnz[i];
1497:       PetscBinaryRead(fd,cols,nz,PETSC_INT);
1498:       MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
1499:     }
1500:     PetscFree(cols);
1501:   } else {
1502:     /* determine buffer space needed for message */
1503:     nz = 0;
1504:     for (i=0; i<m; i++) {
1505:       nz += ourlens[i];
1506:     }
1507:     PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);

1509:     /* receive message of column indices*/
1510:     MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
1511:     MPI_Get_count(&status,MPIU_INT,&maxnz);
1512:     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
1513:   }

1515:   /* loop over local rows, determining number of off diagonal entries */
1516:   PetscMemzero(offlens,m*sizeof(PetscInt));
1517:   jj = 0;
1518:   for (i=0; i<m; i++) {
1519:     for (j=0; j<ourlens[i]; j++) {
1520:       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
1521:       jj++;
1522:     }
1523:   }

1525:   /* create our matrix */
1526:   for (i=0; i<m; i++) {
1527:     ourlens[i] -= offlens[i];
1528:   }
1529:   MatCreate(comm,newmat);
1530:   MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);
1531:   MatSetType(*newmat,type);
1532:   MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);
1533:   A = *newmat;
1534:   for (i=0; i<m; i++) {
1535:     ourlens[i] += offlens[i];
1536:   }

1538:   if (!rank) {
1539:     PetscMalloc(maxnz*sizeof(PetscScalar),&vals);

1541:     /* read in my part of the matrix numerical values  */
1542:     nz = procsnz[0];
1543:     PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1544: 
1545:     /* insert into matrix */
1546:     jj      = rstart;
1547:     smycols = mycols;
1548:     svals   = vals;
1549:     for (i=0; i<m; i++) {
1550:       MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1551:       smycols += ourlens[i];
1552:       svals   += ourlens[i];
1553:       jj++;
1554:     }

1556:     /* read in other processors and ship out */
1557:     for (i=1; i<size; i++) {
1558:       nz   = procsnz[i];
1559:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1560:       MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);
1561:     }
1562:     PetscFree(procsnz);
1563:   } else {
1564:     /* receive numeric values */
1565:     PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);

1567:     /* receive message of values*/
1568:     MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);
1569:     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
1570:     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");

1572:     /* insert into matrix */
1573:     jj      = rstart;
1574:     smycols = mycols;
1575:     svals   = vals;
1576:     for (i=0; i<m; i++) {
1577:       MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1578:       smycols += ourlens[i];
1579:       svals   += ourlens[i];
1580:       jj++;
1581:     }
1582:   }
1583:   PetscFree(ourlens);
1584:   PetscFree(vals);
1585:   PetscFree(mycols);
1586:   PetscFree(rowners);

1588:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1589:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1590:   return(0);
1591: }

1595: PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscTruth *flag)
1596: {
1597:   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
1598:   Mat            a,b;
1599:   PetscTruth     flg;

1603:   a = matA->A;
1604:   b = matB->A;
1605:   MatEqual(a,b,&flg);
1606:   return(0);
1607: }