Actual source code: mpidense.c

petsc-3.5.0 2014-06-30
Report Typos and Errors
  2: /*
  3:    Basic functions for basic parallel dense matrices.
  4: */


  7: #include <../src/mat/impls/dense/mpi/mpidense.h>    /*I   "petscmat.h"  I*/
  8: #include <../src/mat/impls/aij/mpi/mpiaij.h>

 12: /*@

 14:       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
 15:               matrix that represents the operator. For sequential matrices it returns itself.

 17:     Input Parameter:
 18: .      A - the Seq or MPI dense matrix

 20:     Output Parameter:
 21: .      B - the inner matrix

 23:     Level: intermediate

 25: @*/
 26: PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
 27: {
 28:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
 30:   PetscBool      flg;

 33:   PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&flg);
 34:   if (flg) *B = mat->A;
 35:   else *B = A;
 36:   return(0);
 37: }

 41: PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
 42: {
 43:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
 45:   PetscInt       lrow,rstart = A->rmap->rstart,rend = A->rmap->rend;

 48:   if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows");
 49:   lrow = row - rstart;
 50:   MatGetRow(mat->A,lrow,nz,(const PetscInt**)idx,(const PetscScalar**)v);
 51:   return(0);
 52: }

 56: PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
 57: {

 61:   if (idx) {PetscFree(*idx);}
 62:   if (v) {PetscFree(*v);}
 63:   return(0);
 64: }

 68: PetscErrorCode  MatGetDiagonalBlock_MPIDense(Mat A,Mat *a)
 69: {
 70:   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
 72:   PetscInt       m = A->rmap->n,rstart = A->rmap->rstart;
 73:   PetscScalar    *array;
 74:   MPI_Comm       comm;
 75:   Mat            B;

 78:   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only square matrices supported.");

 80:   PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);
 81:   if (!B) {
 82:     PetscObjectGetComm((PetscObject)(mdn->A),&comm);
 83:     MatCreate(comm,&B);
 84:     MatSetSizes(B,m,m,m,m);
 85:     MatSetType(B,((PetscObject)mdn->A)->type_name);
 86:     MatDenseGetArray(mdn->A,&array);
 87:     MatSeqDenseSetPreallocation(B,array+m*rstart);
 88:     MatDenseRestoreArray(mdn->A,&array);
 89:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 90:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 91:     PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);
 92:     *a   = B;
 93:     MatDestroy(&B);
 94:   } else *a = B;
 95:   return(0);
 96: }

100: PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
101: {
102:   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
104:   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
105:   PetscBool      roworiented = A->roworiented;

108:   for (i=0; i<m; i++) {
109:     if (idxm[i] < 0) continue;
110:     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
111:     if (idxm[i] >= rstart && idxm[i] < rend) {
112:       row = idxm[i] - rstart;
113:       if (roworiented) {
114:         MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);
115:       } else {
116:         for (j=0; j<n; j++) {
117:           if (idxn[j] < 0) continue;
118:           if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
119:           MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);
120:         }
121:       }
122:     } else if (!A->donotstash) {
123:       mat->assembled = PETSC_FALSE;
124:       if (roworiented) {
125:         MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n,PETSC_FALSE);
126:       } else {
127:         MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m,PETSC_FALSE);
128:       }
129:     }
130:   }
131:   return(0);
132: }

136: PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
137: {
138:   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
140:   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;

143:   for (i=0; i<m; i++) {
144:     if (idxm[i] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
145:     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
146:     if (idxm[i] >= rstart && idxm[i] < rend) {
147:       row = idxm[i] - rstart;
148:       for (j=0; j<n; j++) {
149:         if (idxn[j] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
150:         if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
151:         MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);
152:       }
153:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
154:   }
155:   return(0);
156: }

160: PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[])
161: {
162:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

166:   MatDenseGetArray(a->A,array);
167:   return(0);
168: }

172: static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
173: {
174:   Mat_MPIDense   *mat  = (Mat_MPIDense*)A->data,*newmatd;
175:   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
177:   PetscInt       i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols;
178:   const PetscInt *irow,*icol;
179:   PetscScalar    *av,*bv,*v = lmat->v;
180:   Mat            newmat;
181:   IS             iscol_local;

184:   ISAllGather(iscol,&iscol_local);
185:   ISGetIndices(isrow,&irow);
186:   ISGetIndices(iscol_local,&icol);
187:   ISGetLocalSize(isrow,&nrows);
188:   ISGetLocalSize(iscol,&ncols);
189:   ISGetSize(iscol,&Ncols); /* global number of columns, size of iscol_local */

191:   /* No parallel redistribution currently supported! Should really check each index set
192:      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
193:      original matrix! */

195:   MatGetLocalSize(A,&nlrows,&nlcols);
196:   MatGetOwnershipRange(A,&rstart,&rend);

198:   /* Check submatrix call */
199:   if (scall == MAT_REUSE_MATRIX) {
200:     /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
201:     /* Really need to test rows and column sizes! */
202:     newmat = *B;
203:   } else {
204:     /* Create and fill new matrix */
205:     MatCreate(PetscObjectComm((PetscObject)A),&newmat);
206:     MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);
207:     MatSetType(newmat,((PetscObject)A)->type_name);
208:     MatMPIDenseSetPreallocation(newmat,NULL);
209:   }

211:   /* Now extract the data pointers and do the copy, column at a time */
212:   newmatd = (Mat_MPIDense*)newmat->data;
213:   bv      = ((Mat_SeqDense*)newmatd->A->data)->v;

215:   for (i=0; i<Ncols; i++) {
216:     av = v + ((Mat_SeqDense*)mat->A->data)->lda*icol[i];
217:     for (j=0; j<nrows; j++) {
218:       *bv++ = av[irow[j] - rstart];
219:     }
220:   }

222:   /* Assemble the matrices so that the correct flags are set */
223:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
224:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);

226:   /* Free work space */
227:   ISRestoreIndices(isrow,&irow);
228:   ISRestoreIndices(iscol_local,&icol);
229:   ISDestroy(&iscol_local);
230:   *B   = newmat;
231:   return(0);
232: }

236: PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[])
237: {
238:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

242:   MatDenseRestoreArray(a->A,array);
243:   return(0);
244: }

248: PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
249: {
250:   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
251:   MPI_Comm       comm;
253:   PetscInt       nstash,reallocs;
254:   InsertMode     addv;

257:   PetscObjectGetComm((PetscObject)mat,&comm);
258:   /* make sure all processors are either in INSERTMODE or ADDMODE */
259:   MPI_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);
260:   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
261:   mat->insertmode = addv; /* in case this processor had no cache */

263:   MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
264:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
265:   PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
266:   return(0);
267: }

271: PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
272: {
273:   Mat_MPIDense   *mdn=(Mat_MPIDense*)mat->data;
275:   PetscInt       i,*row,*col,flg,j,rstart,ncols;
276:   PetscMPIInt    n;
277:   PetscScalar    *val;
278:   InsertMode     addv=mat->insertmode;

281:   /*  wait on receives */
282:   while (1) {
283:     MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
284:     if (!flg) break;

286:     for (i=0; i<n;) {
287:       /* Now identify the consecutive vals belonging to the same row */
288:       for (j=i,rstart=row[j]; j<n; j++) {
289:         if (row[j] != rstart) break;
290:       }
291:       if (j < n) ncols = j-i;
292:       else       ncols = n-i;
293:       /* Now assemble all these values with a single function call */
294:       MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);
295:       i    = j;
296:     }
297:   }
298:   MatStashScatterEnd_Private(&mat->stash);

300:   MatAssemblyBegin(mdn->A,mode);
301:   MatAssemblyEnd(mdn->A,mode);

303:   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
304:     MatSetUpMultiply_MPIDense(mat);
305:   }
306:   return(0);
307: }

311: PetscErrorCode MatZeroEntries_MPIDense(Mat A)
312: {
314:   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;

317:   MatZeroEntries(l->A);
318:   return(0);
319: }

321: /* the code does not do the diagonal entries correctly unless the
322:    matrix is square and the column and row owerships are identical.
323:    This is a BUG. The only way to fix it seems to be to access
324:    mdn->A and mdn->B directly and not through the MatZeroRows()
325:    routine.
326: */
329: PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
330: {
331:   Mat_MPIDense      *l = (Mat_MPIDense*)A->data;
332:   PetscErrorCode    ierr;
333:   PetscInt          i,*owners = A->rmap->range;
334:   PetscInt          *sizes,j,idx,nsends;
335:   PetscInt          nmax,*svalues,*starts,*owner,nrecvs;
336:   PetscInt          *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source;
337:   PetscInt          *lens,*lrows,*values;
338:   PetscMPIInt       n,imdex,rank = l->rank,size = l->size;
339:   MPI_Comm          comm;
340:   MPI_Request       *send_waits,*recv_waits;
341:   MPI_Status        recv_status,*send_status;
342:   PetscBool         found;
343:   const PetscScalar *xx;
344:   PetscScalar       *bb;

347:   PetscObjectGetComm((PetscObject)A,&comm);
348:   if (A->rmap->N != A->cmap->N) SETERRQ(comm,PETSC_ERR_SUP,"Only handles square matrices");
349:   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership");
350:   /*  first count number of contributors to each processor */
351:   PetscCalloc1(2*size,&sizes);
352:   PetscMalloc1(N+1,&owner);  /* see note*/
353:   for (i=0; i<N; i++) {
354:     idx   = rows[i];
355:     found = PETSC_FALSE;
356:     for (j=0; j<size; j++) {
357:       if (idx >= owners[j] && idx < owners[j+1]) {
358:         sizes[2*j]++; sizes[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
359:       }
360:     }
361:     if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
362:   }
363:   nsends = 0;
364:   for (i=0; i<size; i++) nsends += sizes[2*i+1];

366:   /* inform other processors of number of messages and max length*/
367:   PetscMaxSum(comm,sizes,&nmax,&nrecvs);

369:   /* post receives:   */
370:   PetscMalloc1((nrecvs+1)*(nmax+1),&rvalues);
371:   PetscMalloc1((nrecvs+1),&recv_waits);
372:   for (i=0; i<nrecvs; i++) {
373:     MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
374:   }

376:   /* do sends:
377:       1) starts[i] gives the starting index in svalues for stuff going to
378:          the ith processor
379:   */
380:   PetscMalloc1((N+1),&svalues);
381:   PetscMalloc1((nsends+1),&send_waits);
382:   PetscMalloc1((size+1),&starts);

384:   starts[0] = 0;
385:   for (i=1; i<size; i++) starts[i] = starts[i-1] + sizes[2*i-2];
386:   for (i=0; i<N; i++) svalues[starts[owner[i]]++] = rows[i];

388:   starts[0] = 0;
389:   for (i=1; i<size+1; i++) starts[i] = starts[i-1] + sizes[2*i-2];
390:   count = 0;
391:   for (i=0; i<size; i++) {
392:     if (sizes[2*i+1]) {
393:       MPI_Isend(svalues+starts[i],sizes[2*i],MPIU_INT,i,tag,comm,send_waits+count++);
394:     }
395:   }
396:   PetscFree(starts);

398:   base = owners[rank];

400:   /*  wait on receives */
401:   PetscMalloc2(nrecvs,&lens,nrecvs,&source);
402:   count = nrecvs;
403:   slen  = 0;
404:   while (count) {
405:     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
406:     /* unpack receives into our local space */
407:     MPI_Get_count(&recv_status,MPIU_INT,&n);

409:     source[imdex] = recv_status.MPI_SOURCE;
410:     lens[imdex]   = n;
411:     slen += n;
412:     count--;
413:   }
414:   PetscFree(recv_waits);

416:   /* move the data into the send scatter */
417:   PetscMalloc1((slen+1),&lrows);
418:   count = 0;
419:   for (i=0; i<nrecvs; i++) {
420:     values = rvalues + i*nmax;
421:     for (j=0; j<lens[i]; j++) {
422:       lrows[count++] = values[j] - base;
423:     }
424:   }
425:   PetscFree(rvalues);
426:   PetscFree2(lens,source);
427:   PetscFree(owner);
428:   PetscFree(sizes);

430:   /* fix right hand side if needed */
431:   if (x && b) {
432:     VecGetArrayRead(x,&xx);
433:     VecGetArray(b,&bb);
434:     for (i=0; i<slen; i++) {
435:       bb[lrows[i]] = diag*xx[lrows[i]];
436:     }
437:     VecRestoreArrayRead(x,&xx);
438:     VecRestoreArray(b,&bb);
439:   }

441:   /* actually zap the local rows */
442:   MatZeroRows(l->A,slen,lrows,0.0,0,0);
443:   if (diag != 0.0) {
444:     Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data;
445:     PetscInt     m   = ll->lda, i;

447:     for (i=0; i<slen; i++) {
448:       ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag;
449:     }
450:   }
451:   PetscFree(lrows);

453:   /* wait on sends */
454:   if (nsends) {
455:     PetscMalloc1(nsends,&send_status);
456:     MPI_Waitall(nsends,send_waits,send_status);
457:     PetscFree(send_status);
458:   }
459:   PetscFree(send_waits);
460:   PetscFree(svalues);
461:   return(0);
462: }

466: PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
467: {
468:   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;

472:   VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);
473:   VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);
474:   MatMult_SeqDense(mdn->A,mdn->lvec,yy);
475:   return(0);
476: }

480: PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
481: {
482:   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;

486:   VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);
487:   VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);
488:   MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);
489:   return(0);
490: }

494: PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
495: {
496:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
498:   PetscScalar    zero = 0.0;

501:   VecSet(yy,zero);
502:   MatMultTranspose_SeqDense(a->A,xx,a->lvec);
503:   VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
504:   VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
505:   return(0);
506: }

510: PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
511: {
512:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

516:   VecCopy(yy,zz);
517:   MatMultTranspose_SeqDense(a->A,xx,a->lvec);
518:   VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
519:   VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
520:   return(0);
521: }

525: PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
526: {
527:   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
528:   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
530:   PetscInt       len,i,n,m = A->rmap->n,radd;
531:   PetscScalar    *x,zero = 0.0;

534:   VecSet(v,zero);
535:   VecGetArray(v,&x);
536:   VecGetSize(v,&n);
537:   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
538:   len  = PetscMin(a->A->rmap->n,a->A->cmap->n);
539:   radd = A->rmap->rstart*m;
540:   for (i=0; i<len; i++) {
541:     x[i] = aloc->v[radd + i*m + i];
542:   }
543:   VecRestoreArray(v,&x);
544:   return(0);
545: }

549: PetscErrorCode MatDestroy_MPIDense(Mat mat)
550: {
551:   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;

555: #if defined(PETSC_USE_LOG)
556:   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N);
557: #endif
558:   MatStashDestroy_Private(&mat->stash);
559:   MatDestroy(&mdn->A);
560:   VecDestroy(&mdn->lvec);
561:   VecScatterDestroy(&mdn->Mvctx);

563:   PetscFree(mat->data);
564:   PetscObjectChangeTypeName((PetscObject)mat,0);
565:   PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",NULL);
566:   PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",NULL);
567:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",NULL);
568:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",NULL);
569:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",NULL);
570:   return(0);
571: }

575: static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
576: {
577:   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
578:   PetscErrorCode    ierr;
579:   PetscViewerFormat format;
580:   int               fd;
581:   PetscInt          header[4],mmax,N = mat->cmap->N,i,j,m,k;
582:   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
583:   PetscScalar       *work,*v,*vv;
584:   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;

587:   if (mdn->size == 1) {
588:     MatView(mdn->A,viewer);
589:   } else {
590:     PetscViewerBinaryGetDescriptor(viewer,&fd);
591:     MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
592:     MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);

594:     PetscViewerGetFormat(viewer,&format);
595:     if (format == PETSC_VIEWER_NATIVE) {

597:       if (!rank) {
598:         /* store the matrix as a dense matrix */
599:         header[0] = MAT_FILE_CLASSID;
600:         header[1] = mat->rmap->N;
601:         header[2] = N;
602:         header[3] = MATRIX_BINARY_FORMAT_DENSE;
603:         PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);

605:         /* get largest work array needed for transposing array */
606:         mmax = mat->rmap->n;
607:         for (i=1; i<size; i++) {
608:           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
609:         }
610:         PetscMalloc1(mmax*N,&work);

612:         /* write out local array, by rows */
613:         m = mat->rmap->n;
614:         v = a->v;
615:         for (j=0; j<N; j++) {
616:           for (i=0; i<m; i++) {
617:             work[j + i*N] = *v++;
618:           }
619:         }
620:         PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);
621:         /* get largest work array to receive messages from other processes, excludes process zero */
622:         mmax = 0;
623:         for (i=1; i<size; i++) {
624:           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
625:         }
626:         PetscMalloc1(mmax*N,&vv);
627:         for (k = 1; k < size; k++) {
628:           v    = vv;
629:           m    = mat->rmap->range[k+1] - mat->rmap->range[k];
630:           MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,PetscObjectComm((PetscObject)mat));

632:           for (j = 0; j < N; j++) {
633:             for (i = 0; i < m; i++) {
634:               work[j + i*N] = *v++;
635:             }
636:           }
637:           PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);
638:         }
639:         PetscFree(work);
640:         PetscFree(vv);
641:       } else {
642:         MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));
643:       }
644:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE)");
645:   }
646:   return(0);
647: }

649: extern PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
650: #include <petscdraw.h>
653: static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
654: {
655:   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
656:   PetscErrorCode    ierr;
657:   PetscMPIInt       rank = mdn->rank;
658:   PetscViewerType   vtype;
659:   PetscBool         iascii,isdraw;
660:   PetscViewer       sviewer;
661:   PetscViewerFormat format;

664:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
665:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
666:   if (iascii) {
667:     PetscViewerGetType(viewer,&vtype);
668:     PetscViewerGetFormat(viewer,&format);
669:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
670:       MatInfo info;
671:       MatGetInfo(mat,MAT_LOCAL,&info);
672:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
673:       PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);
674:       PetscViewerFlush(viewer);
675:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
676:       VecScatterView(mdn->Mvctx,viewer);
677:       return(0);
678:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
679:       return(0);
680:     }
681:   } else if (isdraw) {
682:     PetscDraw draw;
683:     PetscBool isnull;

685:     PetscViewerDrawGetDraw(viewer,0,&draw);
686:     PetscDrawIsNull(draw,&isnull);
687:     if (isnull) return(0);
688:   }

690:   {
691:     /* assemble the entire matrix onto first processor. */
692:     Mat         A;
693:     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
694:     PetscInt    *cols;
695:     PetscScalar *vals;

697:     MatCreate(PetscObjectComm((PetscObject)mat),&A);
698:     if (!rank) {
699:       MatSetSizes(A,M,N,M,N);
700:     } else {
701:       MatSetSizes(A,0,0,M,N);
702:     }
703:     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
704:     MatSetType(A,MATMPIDENSE);
705:     MatMPIDenseSetPreallocation(A,NULL);
706:     PetscLogObjectParent((PetscObject)mat,(PetscObject)A);

708:     /* Copy the matrix ... This isn't the most efficient means,
709:        but it's quick for now */
710:     A->insertmode = INSERT_VALUES;

712:     row = mat->rmap->rstart;
713:     m   = mdn->A->rmap->n;
714:     for (i=0; i<m; i++) {
715:       MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);
716:       MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);
717:       MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);
718:       row++;
719:     }

721:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
722:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
723:     PetscViewerGetSingleton(viewer,&sviewer);
724:     if (!rank) {
725:         MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);
726:     }
727:     PetscViewerRestoreSingleton(viewer,&sviewer);
728:     PetscViewerFlush(viewer);
729:     MatDestroy(&A);
730:   }
731:   return(0);
732: }

736: PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
737: {
739:   PetscBool      iascii,isbinary,isdraw,issocket;

742:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
743:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
744:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
745:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);

747:   if (iascii || issocket || isdraw) {
748:     MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);
749:   } else if (isbinary) {
750:     MatView_MPIDense_Binary(mat,viewer);
751:   }
752:   return(0);
753: }

757: PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
758: {
759:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
760:   Mat            mdn  = mat->A;
762:   PetscReal      isend[5],irecv[5];

765:   info->block_size = 1.0;

767:   MatGetInfo(mdn,MAT_LOCAL,info);

769:   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
770:   isend[3] = info->memory;  isend[4] = info->mallocs;
771:   if (flag == MAT_LOCAL) {
772:     info->nz_used      = isend[0];
773:     info->nz_allocated = isend[1];
774:     info->nz_unneeded  = isend[2];
775:     info->memory       = isend[3];
776:     info->mallocs      = isend[4];
777:   } else if (flag == MAT_GLOBAL_MAX) {
778:     MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));

780:     info->nz_used      = irecv[0];
781:     info->nz_allocated = irecv[1];
782:     info->nz_unneeded  = irecv[2];
783:     info->memory       = irecv[3];
784:     info->mallocs      = irecv[4];
785:   } else if (flag == MAT_GLOBAL_SUM) {
786:     MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));

788:     info->nz_used      = irecv[0];
789:     info->nz_allocated = irecv[1];
790:     info->nz_unneeded  = irecv[2];
791:     info->memory       = irecv[3];
792:     info->mallocs      = irecv[4];
793:   }
794:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
795:   info->fill_ratio_needed = 0;
796:   info->factor_mallocs    = 0;
797:   return(0);
798: }

802: PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg)
803: {
804:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

808:   switch (op) {
809:   case MAT_NEW_NONZERO_LOCATIONS:
810:   case MAT_NEW_NONZERO_LOCATION_ERR:
811:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
812:     MatSetOption(a->A,op,flg);
813:     break;
814:   case MAT_ROW_ORIENTED:
815:     a->roworiented = flg;

817:     MatSetOption(a->A,op,flg);
818:     break;
819:   case MAT_NEW_DIAGONALS:
820:   case MAT_KEEP_NONZERO_PATTERN:
821:   case MAT_USE_HASH_TABLE:
822:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
823:     break;
824:   case MAT_IGNORE_OFF_PROC_ENTRIES:
825:     a->donotstash = flg;
826:     break;
827:   case MAT_SYMMETRIC:
828:   case MAT_STRUCTURALLY_SYMMETRIC:
829:   case MAT_HERMITIAN:
830:   case MAT_SYMMETRY_ETERNAL:
831:   case MAT_IGNORE_LOWER_TRIANGULAR:
832:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
833:     break;
834:   default:
835:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
836:   }
837:   return(0);
838: }


843: PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
844: {
845:   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
846:   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
847:   PetscScalar    *l,*r,x,*v;
849:   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n;

852:   MatGetLocalSize(A,&s2,&s3);
853:   if (ll) {
854:     VecGetLocalSize(ll,&s2a);
855:     if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
856:     VecGetArray(ll,&l);
857:     for (i=0; i<m; i++) {
858:       x = l[i];
859:       v = mat->v + i;
860:       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
861:     }
862:     VecRestoreArray(ll,&l);
863:     PetscLogFlops(n*m);
864:   }
865:   if (rr) {
866:     VecGetLocalSize(rr,&s3a);
867:     if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
868:     VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);
869:     VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);
870:     VecGetArray(mdn->lvec,&r);
871:     for (i=0; i<n; i++) {
872:       x = r[i];
873:       v = mat->v + i*m;
874:       for (j=0; j<m; j++) (*v++) *= x;
875:     }
876:     VecRestoreArray(mdn->lvec,&r);
877:     PetscLogFlops(n*m);
878:   }
879:   return(0);
880: }

884: PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
885: {
886:   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
887:   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
889:   PetscInt       i,j;
890:   PetscReal      sum = 0.0;
891:   PetscScalar    *v  = mat->v;

894:   if (mdn->size == 1) {
895:      MatNorm(mdn->A,type,nrm);
896:   } else {
897:     if (type == NORM_FROBENIUS) {
898:       for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) {
899:         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
900:       }
901:       MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));
902:       *nrm = PetscSqrtReal(*nrm);
903:       PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);
904:     } else if (type == NORM_1) {
905:       PetscReal *tmp,*tmp2;
906:       PetscMalloc2(A->cmap->N,&tmp,A->cmap->N,&tmp2);
907:       PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));
908:       PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));
909:       *nrm = 0.0;
910:       v    = mat->v;
911:       for (j=0; j<mdn->A->cmap->n; j++) {
912:         for (i=0; i<mdn->A->rmap->n; i++) {
913:           tmp[j] += PetscAbsScalar(*v);  v++;
914:         }
915:       }
916:       MPI_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));
917:       for (j=0; j<A->cmap->N; j++) {
918:         if (tmp2[j] > *nrm) *nrm = tmp2[j];
919:       }
920:       PetscFree2(tmp,tmp);
921:       PetscLogFlops(A->cmap->n*A->rmap->n);
922:     } else if (type == NORM_INFINITY) { /* max row norm */
923:       PetscReal ntemp;
924:       MatNorm(mdn->A,type,&ntemp);
925:       MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));
926:     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm");
927:   }
928:   return(0);
929: }

933: PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
934: {
935:   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
936:   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
937:   Mat            B;
938:   PetscInt       M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart;
940:   PetscInt       j,i;
941:   PetscScalar    *v;

944:   if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place");
945:   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
946:     MatCreate(PetscObjectComm((PetscObject)A),&B);
947:     MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);
948:     MatSetType(B,((PetscObject)A)->type_name);
949:     MatMPIDenseSetPreallocation(B,NULL);
950:   } else {
951:     B = *matout;
952:   }

954:   m    = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v;
955:   PetscMalloc1(m,&rwork);
956:   for (i=0; i<m; i++) rwork[i] = rstart + i;
957:   for (j=0; j<n; j++) {
958:     MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);
959:     v   += m;
960:   }
961:   PetscFree(rwork);
962:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
963:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
964:   if (reuse == MAT_INITIAL_MATRIX || *matout != A) {
965:     *matout = B;
966:   } else {
967:     MatHeaderMerge(A,B);
968:   }
969:   return(0);
970: }


973: static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*);
974: extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar);

978: PetscErrorCode MatSetUp_MPIDense(Mat A)
979: {

983:    MatMPIDenseSetPreallocation(A,0);
984:   return(0);
985: }

989: PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
990: {
992:   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;

995:   MatAXPY(A->A,alpha,B->A,str);
996:   PetscObjectStateIncrease((PetscObject)Y);
997:   return(0);
998: }

1002: PetscErrorCode  MatConjugate_MPIDense(Mat mat)
1003: {
1004:   Mat_MPIDense   *a = (Mat_MPIDense*)mat->data;

1008:   MatConjugate(a->A);
1009:   return(0);
1010: }

1014: PetscErrorCode MatRealPart_MPIDense(Mat A)
1015: {
1016:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1020:   MatRealPart(a->A);
1021:   return(0);
1022: }

1026: PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1027: {
1028:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1032:   MatImaginaryPart(a->A);
1033:   return(0);
1034: }

1036: extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);
1039: PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
1040: {
1042:   PetscInt       i,n;
1043:   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
1044:   PetscReal      *work;

1047:   MatGetSize(A,NULL,&n);
1048:   PetscMalloc1(n,&work);
1049:   MatGetColumnNorms_SeqDense(a->A,type,work);
1050:   if (type == NORM_2) {
1051:     for (i=0; i<n; i++) work[i] *= work[i];
1052:   }
1053:   if (type == NORM_INFINITY) {
1054:     MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);
1055:   } else {
1056:     MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);
1057:   }
1058:   PetscFree(work);
1059:   if (type == NORM_2) {
1060:     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
1061:   }
1062:   return(0);
1063: }

1067: static PetscErrorCode  MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
1068: {
1069:   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;
1071:   PetscScalar    *a;
1072:   PetscInt       m,n,i;

1075:   MatGetSize(d->A,&m,&n);
1076:   MatDenseGetArray(d->A,&a);
1077:   for (i=0; i<m*n; i++) {
1078:     PetscRandomGetValue(rctx,a+i);
1079:   }
1080:   MatDenseRestoreArray(d->A,&a);
1081:   return(0);
1082: }

1084: /* -------------------------------------------------------------------*/
1085: static struct _MatOps MatOps_Values = { MatSetValues_MPIDense,
1086:                                         MatGetRow_MPIDense,
1087:                                         MatRestoreRow_MPIDense,
1088:                                         MatMult_MPIDense,
1089:                                 /*  4*/ MatMultAdd_MPIDense,
1090:                                         MatMultTranspose_MPIDense,
1091:                                         MatMultTransposeAdd_MPIDense,
1092:                                         0,
1093:                                         0,
1094:                                         0,
1095:                                 /* 10*/ 0,
1096:                                         0,
1097:                                         0,
1098:                                         0,
1099:                                         MatTranspose_MPIDense,
1100:                                 /* 15*/ MatGetInfo_MPIDense,
1101:                                         MatEqual_MPIDense,
1102:                                         MatGetDiagonal_MPIDense,
1103:                                         MatDiagonalScale_MPIDense,
1104:                                         MatNorm_MPIDense,
1105:                                 /* 20*/ MatAssemblyBegin_MPIDense,
1106:                                         MatAssemblyEnd_MPIDense,
1107:                                         MatSetOption_MPIDense,
1108:                                         MatZeroEntries_MPIDense,
1109:                                 /* 24*/ MatZeroRows_MPIDense,
1110:                                         0,
1111:                                         0,
1112:                                         0,
1113:                                         0,
1114:                                 /* 29*/ MatSetUp_MPIDense,
1115:                                         0,
1116:                                         0,
1117:                                         0,
1118:                                         0,
1119:                                 /* 34*/ MatDuplicate_MPIDense,
1120:                                         0,
1121:                                         0,
1122:                                         0,
1123:                                         0,
1124:                                 /* 39*/ MatAXPY_MPIDense,
1125:                                         MatGetSubMatrices_MPIDense,
1126:                                         0,
1127:                                         MatGetValues_MPIDense,
1128:                                         0,
1129:                                 /* 44*/ 0,
1130:                                         MatScale_MPIDense,
1131:                                         0,
1132:                                         0,
1133:                                         0,
1134:                                 /* 49*/ MatSetRandom_MPIDense,
1135:                                         0,
1136:                                         0,
1137:                                         0,
1138:                                         0,
1139:                                 /* 54*/ 0,
1140:                                         0,
1141:                                         0,
1142:                                         0,
1143:                                         0,
1144:                                 /* 59*/ MatGetSubMatrix_MPIDense,
1145:                                         MatDestroy_MPIDense,
1146:                                         MatView_MPIDense,
1147:                                         0,
1148:                                         0,
1149:                                 /* 64*/ 0,
1150:                                         0,
1151:                                         0,
1152:                                         0,
1153:                                         0,
1154:                                 /* 69*/ 0,
1155:                                         0,
1156:                                         0,
1157:                                         0,
1158:                                         0,
1159:                                 /* 74*/ 0,
1160:                                         0,
1161:                                         0,
1162:                                         0,
1163:                                         0,
1164:                                 /* 79*/ 0,
1165:                                         0,
1166:                                         0,
1167:                                         0,
1168:                                 /* 83*/ MatLoad_MPIDense,
1169:                                         0,
1170:                                         0,
1171:                                         0,
1172:                                         0,
1173:                                         0,
1174:                                 /* 89*/
1175:                                         0,
1176:                                         0,
1177:                                         0,
1178:                                         0,
1179:                                         0,
1180:                                 /* 94*/ 0,
1181:                                         0,
1182:                                         0,
1183:                                         0,
1184:                                         0,
1185:                                 /* 99*/ 0,
1186:                                         0,
1187:                                         0,
1188:                                         MatConjugate_MPIDense,
1189:                                         0,
1190:                                 /*104*/ 0,
1191:                                         MatRealPart_MPIDense,
1192:                                         MatImaginaryPart_MPIDense,
1193:                                         0,
1194:                                         0,
1195:                                 /*109*/ 0,
1196:                                         0,
1197:                                         0,
1198:                                         0,
1199:                                         0,
1200:                                 /*114*/ 0,
1201:                                         0,
1202:                                         0,
1203:                                         0,
1204:                                         0,
1205:                                 /*119*/ 0,
1206:                                         0,
1207:                                         0,
1208:                                         0,
1209:                                         0,
1210:                                 /*124*/ 0,
1211:                                         MatGetColumnNorms_MPIDense,
1212:                                         0,
1213:                                         0,
1214:                                         0,
1215:                                 /*129*/ 0,
1216:                                         0,
1217:                                         0,
1218:                                         0,
1219:                                         0,
1220:                                 /*134*/ 0,
1221:                                         0,
1222:                                         0,
1223:                                         0,
1224:                                         0,
1225:                                 /*139*/ 0,
1226:                                         0,
1227:                                         0
1228: };

1232: PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1233: {
1234:   Mat_MPIDense   *a;

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

1242:   a       = (Mat_MPIDense*)mat->data;
1243:   PetscLayoutSetUp(mat->rmap);
1244:   PetscLayoutSetUp(mat->cmap);
1245:   a->nvec = mat->cmap->n;

1247:   MatCreate(PETSC_COMM_SELF,&a->A);
1248:   MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);
1249:   MatSetType(a->A,MATSEQDENSE);
1250:   MatSeqDenseSetPreallocation(a->A,data);
1251:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
1252:   return(0);
1253: }

1257: PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat)
1258: {
1259:   Mat_MPIDense   *a;

1263:   PetscNewLog(mat,&a);
1264:   mat->data = (void*)a;
1265:   PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));

1267:   mat->insertmode = NOT_SET_VALUES;
1268:   MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);
1269:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);

1271:   /* build cache for off array entries formed */
1272:   a->donotstash = PETSC_FALSE;

1274:   MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);

1276:   /* stuff used for matrix vector multiply */
1277:   a->lvec        = 0;
1278:   a->Mvctx       = 0;
1279:   a->roworiented = PETSC_TRUE;

1281:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);
1282:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);

1284:   PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPIDense);
1285:   PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);
1286:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",MatMatMult_MPIAIJ_MPIDense);
1287:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",MatMatMultSymbolic_MPIAIJ_MPIDense);
1288:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",MatMatMultNumeric_MPIAIJ_MPIDense);

1290:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",MatTransposeMatMult_MPIAIJ_MPIDense);
1291:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",MatTransposeMatMultSymbolic_MPIAIJ_MPIDense);
1292:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",MatTransposeMatMultNumeric_MPIAIJ_MPIDense);
1293:   PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);
1294:   return(0);
1295: }

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

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

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

1306:   Level: beginner


1309: .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1310: M*/

1314: /*@C
1315:    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries

1317:    Not collective

1319:    Input Parameters:
1320: .  A - the matrix
1321: -  data - optional location of matrix data.  Set data=NULL for PETSc
1322:    to control all matrix memory allocation.

1324:    Notes:
1325:    The dense format is fully compatible with standard Fortran 77
1326:    storage by columns.

1328:    The data input variable is intended primarily for Fortran programmers
1329:    who wish to allocate their own matrix memory space.  Most users should
1330:    set data=NULL.

1332:    Level: intermediate

1334: .keywords: matrix,dense, parallel

1336: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1337: @*/
1338: PetscErrorCode  MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1339: {

1343:   PetscTryMethod(mat,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(mat,data));
1344:   return(0);
1345: }

1349: /*@C
1350:    MatCreateDense - Creates a parallel matrix in dense format.

1352:    Collective on MPI_Comm

1354:    Input Parameters:
1355: +  comm - MPI communicator
1356: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1357: .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
1358: .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1359: .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1360: -  data - optional location of matrix data.  Set data=NULL (NULL_SCALAR for Fortran users) for PETSc
1361:    to control all matrix memory allocation.

1363:    Output Parameter:
1364: .  A - the matrix

1366:    Notes:
1367:    The dense format is fully compatible with standard Fortran 77
1368:    storage by columns.

1370:    The data input variable is intended primarily for Fortran programmers
1371:    who wish to allocate their own matrix memory space.  Most users should
1372:    set data=NULL (NULL_SCALAR for Fortran users).

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

1377:    Level: intermediate

1379: .keywords: matrix,dense, parallel

1381: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1382: @*/
1383: PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
1384: {
1386:   PetscMPIInt    size;

1389:   MatCreate(comm,A);
1390:   MatSetSizes(*A,m,n,M,N);
1391:   MPI_Comm_size(comm,&size);
1392:   if (size > 1) {
1393:     MatSetType(*A,MATMPIDENSE);
1394:     MatMPIDenseSetPreallocation(*A,data);
1395:   } else {
1396:     MatSetType(*A,MATSEQDENSE);
1397:     MatSeqDenseSetPreallocation(*A,data);
1398:   }
1399:   return(0);
1400: }

1404: static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
1405: {
1406:   Mat            mat;
1407:   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;

1411:   *newmat = 0;
1412:   MatCreate(PetscObjectComm((PetscObject)A),&mat);
1413:   MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
1414:   MatSetType(mat,((PetscObject)A)->type_name);
1415:   a       = (Mat_MPIDense*)mat->data;
1416:   PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));

1418:   mat->factortype   = A->factortype;
1419:   mat->assembled    = PETSC_TRUE;
1420:   mat->preallocated = PETSC_TRUE;

1422:   a->size         = oldmat->size;
1423:   a->rank         = oldmat->rank;
1424:   mat->insertmode = NOT_SET_VALUES;
1425:   a->nvec         = oldmat->nvec;
1426:   a->donotstash   = oldmat->donotstash;

1428:   PetscLayoutReference(A->rmap,&mat->rmap);
1429:   PetscLayoutReference(A->cmap,&mat->cmap);

1431:   MatSetUpMultiply_MPIDense(mat);
1432:   MatDuplicate(oldmat->A,cpvalues,&a->A);
1433:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);

1435:   *newmat = mat;
1436:   return(0);
1437: }

1441: PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat,PetscInt sizesset)
1442: {
1444:   PetscMPIInt    rank,size;
1445:   PetscInt       *rowners,i,m,nz,j;
1446:   PetscScalar    *array,*vals,*vals_ptr;

1449:   MPI_Comm_rank(comm,&rank);
1450:   MPI_Comm_size(comm,&size);

1452:   /* determine ownership of all rows */
1453:   if (newmat->rmap->n < 0) m = M/size + ((M % size) > rank);
1454:   else m = newmat->rmap->n;
1455:   PetscMalloc1((size+2),&rowners);
1456:   MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);
1457:   rowners[0] = 0;
1458:   for (i=2; i<=size; i++) {
1459:     rowners[i] += rowners[i-1];
1460:   }

1462:   if (!sizesset) {
1463:     MatSetSizes(newmat,m,PETSC_DECIDE,M,N);
1464:   }
1465:   MatMPIDenseSetPreallocation(newmat,NULL);
1466:   MatDenseGetArray(newmat,&array);

1468:   if (!rank) {
1469:     PetscMalloc1(m*N,&vals);

1471:     /* read in my part of the matrix numerical values  */
1472:     PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);

1474:     /* insert into matrix-by row (this is why cannot directly read into array */
1475:     vals_ptr = vals;
1476:     for (i=0; i<m; i++) {
1477:       for (j=0; j<N; j++) {
1478:         array[i + j*m] = *vals_ptr++;
1479:       }
1480:     }

1482:     /* read in other processors and ship out */
1483:     for (i=1; i<size; i++) {
1484:       nz   = (rowners[i+1] - rowners[i])*N;
1485:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1486:       MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);
1487:     }
1488:   } else {
1489:     /* receive numeric values */
1490:     PetscMalloc1(m*N,&vals);

1492:     /* receive message of values*/
1493:     MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);

1495:     /* insert into matrix-by row (this is why cannot directly read into array */
1496:     vals_ptr = vals;
1497:     for (i=0; i<m; i++) {
1498:       for (j=0; j<N; j++) {
1499:         array[i + j*m] = *vals_ptr++;
1500:       }
1501:     }
1502:   }
1503:   MatDenseRestoreArray(newmat,&array);
1504:   PetscFree(rowners);
1505:   PetscFree(vals);
1506:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1507:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1508:   return(0);
1509: }

1513: PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer)
1514: {
1515:   PetscScalar    *vals,*svals;
1516:   MPI_Comm       comm;
1517:   MPI_Status     status;
1518:   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
1519:   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
1520:   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
1521:   PetscInt       i,nz,j,rstart,rend,sizesset=1,grows,gcols;
1522:   int            fd;

1526:   PetscObjectGetComm((PetscObject)viewer,&comm);
1527:   MPI_Comm_size(comm,&size);
1528:   MPI_Comm_rank(comm,&rank);
1529:   if (!rank) {
1530:     PetscViewerBinaryGetDescriptor(viewer,&fd);
1531:     PetscBinaryRead(fd,(char*)header,4,PETSC_INT);
1532:     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
1533:   }
1534:   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) sizesset = 0;

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

1539:   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
1540:   if (sizesset && newmat->rmap->N < 0) newmat->rmap->N = M;
1541:   if (sizesset && newmat->cmap->N < 0) newmat->cmap->N = N;

1543:   /* If global sizes are set, check if they are consistent with that given in the file */
1544:   if (sizesset) {
1545:     MatGetSize(newmat,&grows,&gcols);
1546:   }
1547:   if (sizesset && newmat->rmap->N != grows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%d) and input matrix has (%d)",M,grows);
1548:   if (sizesset && newmat->cmap->N != gcols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%d) and input matrix has (%d)",N,gcols);

1550:   /*
1551:        Handle case where matrix is stored on disk as a dense matrix
1552:   */
1553:   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1554:     MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat,sizesset);
1555:     return(0);
1556:   }

1558:   /* determine ownership of all rows */
1559:   if (newmat->rmap->n < 0) {
1560:     PetscMPIIntCast(M/size + ((M % size) > rank),&m);
1561:   } else {
1562:     PetscMPIIntCast(newmat->rmap->n,&m);
1563:   }
1564:   PetscMalloc1((size+2),&rowners);
1565:   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1566:   rowners[0] = 0;
1567:   for (i=2; i<=size; i++) {
1568:     rowners[i] += rowners[i-1];
1569:   }
1570:   rstart = rowners[rank];
1571:   rend   = rowners[rank+1];

1573:   /* distribute row lengths to all processors */
1574:   PetscMalloc2(rend-rstart,&ourlens,rend-rstart,&offlens);
1575:   if (!rank) {
1576:     PetscMalloc1(M,&rowlengths);
1577:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
1578:     PetscMalloc1(size,&sndcounts);
1579:     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1580:     MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1581:     PetscFree(sndcounts);
1582:   } else {
1583:     MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1584:   }

1586:   if (!rank) {
1587:     /* calculate the number of nonzeros on each processor */
1588:     PetscMalloc1(size,&procsnz);
1589:     PetscMemzero(procsnz,size*sizeof(PetscInt));
1590:     for (i=0; i<size; i++) {
1591:       for (j=rowners[i]; j< rowners[i+1]; j++) {
1592:         procsnz[i] += rowlengths[j];
1593:       }
1594:     }
1595:     PetscFree(rowlengths);

1597:     /* determine max buffer needed and allocate it */
1598:     maxnz = 0;
1599:     for (i=0; i<size; i++) {
1600:       maxnz = PetscMax(maxnz,procsnz[i]);
1601:     }
1602:     PetscMalloc1(maxnz,&cols);

1604:     /* read in my part of the matrix column indices  */
1605:     nz   = procsnz[0];
1606:     PetscMalloc1(nz,&mycols);
1607:     PetscBinaryRead(fd,mycols,nz,PETSC_INT);

1609:     /* read in every one elses and ship off */
1610:     for (i=1; i<size; i++) {
1611:       nz   = procsnz[i];
1612:       PetscBinaryRead(fd,cols,nz,PETSC_INT);
1613:       MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
1614:     }
1615:     PetscFree(cols);
1616:   } else {
1617:     /* determine buffer space needed for message */
1618:     nz = 0;
1619:     for (i=0; i<m; i++) {
1620:       nz += ourlens[i];
1621:     }
1622:     PetscMalloc1((nz+1),&mycols);

1624:     /* receive message of column indices*/
1625:     MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
1626:     MPI_Get_count(&status,MPIU_INT,&maxnz);
1627:     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
1628:   }

1630:   /* loop over local rows, determining number of off diagonal entries */
1631:   PetscMemzero(offlens,m*sizeof(PetscInt));
1632:   jj   = 0;
1633:   for (i=0; i<m; i++) {
1634:     for (j=0; j<ourlens[i]; j++) {
1635:       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
1636:       jj++;
1637:     }
1638:   }

1640:   /* create our matrix */
1641:   for (i=0; i<m; i++) ourlens[i] -= offlens[i];

1643:   if (!sizesset) {
1644:     MatSetSizes(newmat,m,PETSC_DECIDE,M,N);
1645:   }
1646:   MatMPIDenseSetPreallocation(newmat,NULL);
1647:   for (i=0; i<m; i++) ourlens[i] += offlens[i];

1649:   if (!rank) {
1650:     PetscMalloc1(maxnz,&vals);

1652:     /* read in my part of the matrix numerical values  */
1653:     nz   = procsnz[0];
1654:     PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);

1656:     /* insert into matrix */
1657:     jj      = rstart;
1658:     smycols = mycols;
1659:     svals   = vals;
1660:     for (i=0; i<m; i++) {
1661:       MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1662:       smycols += ourlens[i];
1663:       svals   += ourlens[i];
1664:       jj++;
1665:     }

1667:     /* read in other processors and ship out */
1668:     for (i=1; i<size; i++) {
1669:       nz   = procsnz[i];
1670:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1671:       MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);
1672:     }
1673:     PetscFree(procsnz);
1674:   } else {
1675:     /* receive numeric values */
1676:     PetscMalloc1((nz+1),&vals);

1678:     /* receive message of values*/
1679:     MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);
1680:     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
1681:     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");

1683:     /* insert into matrix */
1684:     jj      = rstart;
1685:     smycols = mycols;
1686:     svals   = vals;
1687:     for (i=0; i<m; i++) {
1688:       MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1689:       smycols += ourlens[i];
1690:       svals   += ourlens[i];
1691:       jj++;
1692:     }
1693:   }
1694:   PetscFree2(ourlens,offlens);
1695:   PetscFree(vals);
1696:   PetscFree(mycols);
1697:   PetscFree(rowners);

1699:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1700:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1701:   return(0);
1702: }

1706: PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
1707: {
1708:   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
1709:   Mat            a,b;
1710:   PetscBool      flg;

1714:   a    = matA->A;
1715:   b    = matB->A;
1716:   MatEqual(a,b,&flg);
1717:   MPI_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1718:   return(0);
1719: }