Actual source code: mpidense.c

petsc-master 2016-02-10
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>
  9: #include <petscblaslapack.h>

 13: /*@

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

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

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

 24:     Level: intermediate

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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,mat->insertmode);
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);

566:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);
567:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);
568: #if defined(PETSC_HAVE_ELEMENTAL)
569:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",NULL);
570: #endif
571:   PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",NULL);
572:   PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",NULL);
573:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",NULL);
574:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",NULL);
575:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",NULL);
576:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",NULL);
577:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",NULL);
578:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",NULL);
579:   return(0);
580: }

584: static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
585: {
586:   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
587:   PetscErrorCode    ierr;
588:   PetscViewerFormat format;
589:   int               fd;
590:   PetscInt          header[4],mmax,N = mat->cmap->N,i,j,m,k;
591:   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
592:   PetscScalar       *work,*v,*vv;
593:   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;

596:   if (mdn->size == 1) {
597:     MatView(mdn->A,viewer);
598:   } else {
599:     PetscViewerBinaryGetDescriptor(viewer,&fd);
600:     MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
601:     MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);

603:     PetscViewerGetFormat(viewer,&format);
604:     if (format == PETSC_VIEWER_NATIVE) {

606:       if (!rank) {
607:         /* store the matrix as a dense matrix */
608:         header[0] = MAT_FILE_CLASSID;
609:         header[1] = mat->rmap->N;
610:         header[2] = N;
611:         header[3] = MATRIX_BINARY_FORMAT_DENSE;
612:         PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);

614:         /* get largest work array needed for transposing array */
615:         mmax = mat->rmap->n;
616:         for (i=1; i<size; i++) {
617:           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
618:         }
619:         PetscMalloc1(mmax*N,&work);

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

641:           for (j = 0; j < N; j++) {
642:             for (i = 0; i < m; i++) {
643:               work[j + i*N] = *v++;
644:             }
645:           }
646:           PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);
647:         }
648:         PetscFree(work);
649:         PetscFree(vv);
650:       } else {
651:         MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));
652:       }
653:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE)");
654:   }
655:   return(0);
656: }

658: extern PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
659: #include <petscdraw.h>
662: static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
663: {
664:   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
665:   PetscErrorCode    ierr;
666:   PetscMPIInt       rank = mdn->rank;
667:   PetscViewerType   vtype;
668:   PetscBool         iascii,isdraw;
669:   PetscViewer       sviewer;
670:   PetscViewerFormat format;

673:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
674:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
675:   if (iascii) {
676:     PetscViewerGetType(viewer,&vtype);
677:     PetscViewerGetFormat(viewer,&format);
678:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
679:       MatInfo info;
680:       MatGetInfo(mat,MAT_LOCAL,&info);
681:       PetscViewerASCIIPushSynchronized(viewer);
682:       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);
683:       PetscViewerFlush(viewer);
684:       PetscViewerASCIIPopSynchronized(viewer);
685:       VecScatterView(mdn->Mvctx,viewer);
686:       return(0);
687:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
688:       return(0);
689:     }
690:   } else if (isdraw) {
691:     PetscDraw draw;
692:     PetscBool isnull;

694:     PetscViewerDrawGetDraw(viewer,0,&draw);
695:     PetscDrawIsNull(draw,&isnull);
696:     if (isnull) return(0);
697:   }

699:   {
700:     /* assemble the entire matrix onto first processor. */
701:     Mat         A;
702:     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
703:     PetscInt    *cols;
704:     PetscScalar *vals;

706:     MatCreate(PetscObjectComm((PetscObject)mat),&A);
707:     if (!rank) {
708:       MatSetSizes(A,M,N,M,N);
709:     } else {
710:       MatSetSizes(A,0,0,M,N);
711:     }
712:     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
713:     MatSetType(A,MATMPIDENSE);
714:     MatMPIDenseSetPreallocation(A,NULL);
715:     PetscLogObjectParent((PetscObject)mat,(PetscObject)A);

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

721:     row = mat->rmap->rstart;
722:     m   = mdn->A->rmap->n;
723:     for (i=0; i<m; i++) {
724:       MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);
725:       MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);
726:       MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);
727:       row++;
728:     }

730:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
731:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
732:     PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
733:     if (!rank) {
734:         MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);
735:     }
736:     PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
737:     PetscViewerFlush(viewer);
738:     MatDestroy(&A);
739:   }
740:   return(0);
741: }

745: PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
746: {
748:   PetscBool      iascii,isbinary,isdraw,issocket;

751:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
752:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
753:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
754:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);

756:   if (iascii || issocket || isdraw) {
757:     MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);
758:   } else if (isbinary) {
759:     MatView_MPIDense_Binary(mat,viewer);
760:   }
761:   return(0);
762: }

766: PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
767: {
768:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
769:   Mat            mdn  = mat->A;
771:   PetscReal      isend[5],irecv[5];

774:   info->block_size = 1.0;

776:   MatGetInfo(mdn,MAT_LOCAL,info);

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

789:     info->nz_used      = irecv[0];
790:     info->nz_allocated = irecv[1];
791:     info->nz_unneeded  = irecv[2];
792:     info->memory       = irecv[3];
793:     info->mallocs      = irecv[4];
794:   } else if (flag == MAT_GLOBAL_SUM) {
795:     MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));

797:     info->nz_used      = irecv[0];
798:     info->nz_allocated = irecv[1];
799:     info->nz_unneeded  = irecv[2];
800:     info->memory       = irecv[3];
801:     info->mallocs      = irecv[4];
802:   }
803:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
804:   info->fill_ratio_needed = 0;
805:   info->factor_mallocs    = 0;
806:   return(0);
807: }

811: PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg)
812: {
813:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

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


851: PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
852: {
853:   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
854:   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
855:   PetscScalar    *l,*r,x,*v;
857:   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n;

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

892: PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
893: {
894:   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
895:   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
897:   PetscInt       i,j;
898:   PetscReal      sum = 0.0;
899:   PetscScalar    *v  = mat->v;

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

941: PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
942: {
943:   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
944:   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
945:   Mat            B;
946:   PetscInt       M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart;
948:   PetscInt       j,i;
949:   PetscScalar    *v;

952:   if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place");
953:   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
954:     MatCreate(PetscObjectComm((PetscObject)A),&B);
955:     MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);
956:     MatSetType(B,((PetscObject)A)->type_name);
957:     MatMPIDenseSetPreallocation(B,NULL);
958:   } else {
959:     B = *matout;
960:   }

962:   m    = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v;
963:   PetscMalloc1(m,&rwork);
964:   for (i=0; i<m; i++) rwork[i] = rstart + i;
965:   for (j=0; j<n; j++) {
966:     MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);
967:     v   += m;
968:   }
969:   PetscFree(rwork);
970:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
971:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
972:   if (reuse == MAT_INITIAL_MATRIX || *matout != A) {
973:     *matout = B;
974:   } else {
975:     MatHeaderMerge(A,&B);
976:   }
977:   return(0);
978: }


981: static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*);
982: extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar);

986: PetscErrorCode MatSetUp_MPIDense(Mat A)
987: {

991:    MatMPIDenseSetPreallocation(A,0);
992:   return(0);
993: }

997: PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
998: {
1000:   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;

1003:   MatAXPY(A->A,alpha,B->A,str);
1004:   PetscObjectStateIncrease((PetscObject)Y);
1005:   return(0);
1006: }

1010: PetscErrorCode  MatConjugate_MPIDense(Mat mat)
1011: {
1012:   Mat_MPIDense   *a = (Mat_MPIDense*)mat->data;

1016:   MatConjugate(a->A);
1017:   return(0);
1018: }

1022: PetscErrorCode MatRealPart_MPIDense(Mat A)
1023: {
1024:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1028:   MatRealPart(a->A);
1029:   return(0);
1030: }

1034: PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1035: {
1036:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1040:   MatImaginaryPart(a->A);
1041:   return(0);
1042: }

1044: extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);
1047: PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
1048: {
1050:   PetscInt       i,n;
1051:   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
1052:   PetscReal      *work;

1055:   MatGetSize(A,NULL,&n);
1056:   PetscMalloc1(n,&work);
1057:   MatGetColumnNorms_SeqDense(a->A,type,work);
1058:   if (type == NORM_2) {
1059:     for (i=0; i<n; i++) work[i] *= work[i];
1060:   }
1061:   if (type == NORM_INFINITY) {
1062:     MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);
1063:   } else {
1064:     MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);
1065:   }
1066:   PetscFree(work);
1067:   if (type == NORM_2) {
1068:     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
1069:   }
1070:   return(0);
1071: }

1075: static PetscErrorCode  MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
1076: {
1077:   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;
1079:   PetscScalar    *a;
1080:   PetscInt       m,n,i;

1083:   MatGetSize(d->A,&m,&n);
1084:   MatDenseGetArray(d->A,&a);
1085:   for (i=0; i<m*n; i++) {
1086:     PetscRandomGetValue(rctx,a+i);
1087:   }
1088:   MatDenseRestoreArray(d->A,&a);
1089:   return(0);
1090: }

1092: extern PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat,Mat);

1096: static PetscErrorCode MatMissingDiagonal_MPIDense(Mat A,PetscBool  *missing,PetscInt *d)
1097: {
1099:   *missing = PETSC_FALSE;
1100:   return(0);
1101: }

1103: /* -------------------------------------------------------------------*/
1104: static struct _MatOps MatOps_Values = { MatSetValues_MPIDense,
1105:                                         MatGetRow_MPIDense,
1106:                                         MatRestoreRow_MPIDense,
1107:                                         MatMult_MPIDense,
1108:                                 /*  4*/ MatMultAdd_MPIDense,
1109:                                         MatMultTranspose_MPIDense,
1110:                                         MatMultTransposeAdd_MPIDense,
1111:                                         0,
1112:                                         0,
1113:                                         0,
1114:                                 /* 10*/ 0,
1115:                                         0,
1116:                                         0,
1117:                                         0,
1118:                                         MatTranspose_MPIDense,
1119:                                 /* 15*/ MatGetInfo_MPIDense,
1120:                                         MatEqual_MPIDense,
1121:                                         MatGetDiagonal_MPIDense,
1122:                                         MatDiagonalScale_MPIDense,
1123:                                         MatNorm_MPIDense,
1124:                                 /* 20*/ MatAssemblyBegin_MPIDense,
1125:                                         MatAssemblyEnd_MPIDense,
1126:                                         MatSetOption_MPIDense,
1127:                                         MatZeroEntries_MPIDense,
1128:                                 /* 24*/ MatZeroRows_MPIDense,
1129:                                         0,
1130:                                         0,
1131:                                         0,
1132:                                         0,
1133:                                 /* 29*/ MatSetUp_MPIDense,
1134:                                         0,
1135:                                         0,
1136:                                         0,
1137:                                         0,
1138:                                 /* 34*/ MatDuplicate_MPIDense,
1139:                                         0,
1140:                                         0,
1141:                                         0,
1142:                                         0,
1143:                                 /* 39*/ MatAXPY_MPIDense,
1144:                                         MatGetSubMatrices_MPIDense,
1145:                                         0,
1146:                                         MatGetValues_MPIDense,
1147:                                         0,
1148:                                 /* 44*/ 0,
1149:                                         MatScale_MPIDense,
1150:                                         MatShift_Basic,
1151:                                         0,
1152:                                         0,
1153:                                 /* 49*/ MatSetRandom_MPIDense,
1154:                                         0,
1155:                                         0,
1156:                                         0,
1157:                                         0,
1158:                                 /* 54*/ 0,
1159:                                         0,
1160:                                         0,
1161:                                         0,
1162:                                         0,
1163:                                 /* 59*/ MatGetSubMatrix_MPIDense,
1164:                                         MatDestroy_MPIDense,
1165:                                         MatView_MPIDense,
1166:                                         0,
1167:                                         0,
1168:                                 /* 64*/ 0,
1169:                                         0,
1170:                                         0,
1171:                                         0,
1172:                                         0,
1173:                                 /* 69*/ 0,
1174:                                         0,
1175:                                         0,
1176:                                         0,
1177:                                         0,
1178:                                 /* 74*/ 0,
1179:                                         0,
1180:                                         0,
1181:                                         0,
1182:                                         0,
1183:                                 /* 79*/ 0,
1184:                                         0,
1185:                                         0,
1186:                                         0,
1187:                                 /* 83*/ MatLoad_MPIDense,
1188:                                         0,
1189:                                         0,
1190:                                         0,
1191:                                         0,
1192:                                         0,
1193: #if defined(PETSC_HAVE_ELEMENTAL)
1194:                                 /* 89*/ MatMatMult_MPIDense_MPIDense,
1195:                                         MatMatMultSymbolic_MPIDense_MPIDense,
1196: #else
1197:                                 /* 89*/ 0,
1198:                                         0,
1199: #endif
1200:                                         MatMatMultNumeric_MPIDense,
1201:                                         0,
1202:                                         0,
1203:                                 /* 94*/ 0,
1204:                                         0,
1205:                                         0,
1206:                                         0,
1207:                                         0,
1208:                                 /* 99*/ 0,
1209:                                         0,
1210:                                         0,
1211:                                         MatConjugate_MPIDense,
1212:                                         0,
1213:                                 /*104*/ 0,
1214:                                         MatRealPart_MPIDense,
1215:                                         MatImaginaryPart_MPIDense,
1216:                                         0,
1217:                                         0,
1218:                                 /*109*/ 0,
1219:                                         0,
1220:                                         0,
1221:                                         0,
1222:                                         MatMissingDiagonal_MPIDense,
1223:                                 /*114*/ 0,
1224:                                         0,
1225:                                         0,
1226:                                         0,
1227:                                         0,
1228:                                 /*119*/ 0,
1229:                                         0,
1230:                                         0,
1231:                                         0,
1232:                                         0,
1233:                                 /*124*/ 0,
1234:                                         MatGetColumnNorms_MPIDense,
1235:                                         0,
1236:                                         0,
1237:                                         0,
1238:                                 /*129*/ 0,
1239:                                         MatTransposeMatMult_MPIDense_MPIDense,
1240:                                         MatTransposeMatMultSymbolic_MPIDense_MPIDense,
1241:                                         MatTransposeMatMultNumeric_MPIDense_MPIDense,
1242:                                         0,
1243:                                 /*134*/ 0,
1244:                                         0,
1245:                                         0,
1246:                                         0,
1247:                                         0,
1248:                                 /*139*/ 0,
1249:                                         0,
1250:                                         0
1251: };

1255: PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1256: {
1257:   Mat_MPIDense   *a;

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

1265:   a       = (Mat_MPIDense*)mat->data;
1266:   PetscLayoutSetUp(mat->rmap);
1267:   PetscLayoutSetUp(mat->cmap);
1268:   a->nvec = mat->cmap->n;

1270:   MatCreate(PETSC_COMM_SELF,&a->A);
1271:   MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);
1272:   MatSetType(a->A,MATSEQDENSE);
1273:   MatSeqDenseSetPreallocation(a->A,data);
1274:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
1275:   return(0);
1276: }

1278: #if defined(PETSC_HAVE_ELEMENTAL)
1281: PETSC_EXTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
1282: {
1283:   Mat            mat_elemental;
1285:   PetscScalar    *v;
1286:   PetscInt       m=A->rmap->n,N=A->cmap->N,rstart=A->rmap->rstart,i,*rows,*cols;
1287: 
1289:   if (reuse == MAT_REUSE_MATRIX) {
1290:     mat_elemental = *newmat;
1291:     MatZeroEntries(*newmat);
1292:   } else {
1293:     MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
1294:     MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);
1295:     MatSetType(mat_elemental,MATELEMENTAL);
1296:     MatSetUp(mat_elemental);
1297:     MatSetOption(mat_elemental,MAT_ROW_ORIENTED,PETSC_FALSE);
1298:   }

1300:   PetscMalloc2(m,&rows,N,&cols);
1301:   for (i=0; i<N; i++) cols[i] = i;
1302:   for (i=0; i<m; i++) rows[i] = rstart + i;
1303: 
1304:   /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
1305:   MatDenseGetArray(A,&v);
1306:   MatSetValues(mat_elemental,m,rows,N,cols,v,ADD_VALUES);
1307:   MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
1308:   MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
1309:   MatDenseRestoreArray(A,&v);
1310:   PetscFree2(rows,cols);

1312:   if (reuse == MAT_INPLACE_MATRIX) {
1313:     MatHeaderReplace(A,&mat_elemental);
1314:   } else {
1315:     *newmat = mat_elemental;
1316:   }
1317:   return(0);
1318: }
1319: #endif

1323: PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat)
1324: {
1325:   Mat_MPIDense   *a;

1329:   PetscNewLog(mat,&a);
1330:   mat->data = (void*)a;
1331:   PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));

1333:   mat->insertmode = NOT_SET_VALUES;
1334:   MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);
1335:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);

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

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

1342:   /* stuff used for matrix vector multiply */
1343:   a->lvec        = 0;
1344:   a->Mvctx       = 0;
1345:   a->roworiented = PETSC_TRUE;

1347:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);
1348:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);
1349: #if defined(PETSC_HAVE_ELEMENTAL)
1350:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",MatConvert_MPIDense_Elemental);
1351: #endif
1352:   PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPIDense);
1353:   PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);
1354:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",MatMatMult_MPIAIJ_MPIDense);
1355:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",MatMatMultSymbolic_MPIAIJ_MPIDense);
1356:   PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",MatMatMultNumeric_MPIAIJ_MPIDense);

1358:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",MatTransposeMatMult_MPIAIJ_MPIDense);
1359:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",MatTransposeMatMultSymbolic_MPIAIJ_MPIDense);
1360:   PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",MatTransposeMatMultNumeric_MPIAIJ_MPIDense);
1361:   PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);
1362:   return(0);
1363: }

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

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

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

1374:   Level: beginner


1377: .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1378: M*/

1382: /*@C
1383:    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries

1385:    Not collective

1387:    Input Parameters:
1388: .  B - the matrix
1389: -  data - optional location of matrix data.  Set data=NULL for PETSc
1390:    to control all matrix memory allocation.

1392:    Notes:
1393:    The dense format is fully compatible with standard Fortran 77
1394:    storage by columns.

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

1400:    Level: intermediate

1402: .keywords: matrix,dense, parallel

1404: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1405: @*/
1406: PetscErrorCode  MatMPIDenseSetPreallocation(Mat B,PetscScalar *data)
1407: {

1411:   PetscTryMethod(B,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(B,data));
1412:   return(0);
1413: }

1417: /*@C
1418:    MatCreateDense - Creates a parallel matrix in dense format.

1420:    Collective on MPI_Comm

1422:    Input Parameters:
1423: +  comm - MPI communicator
1424: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1425: .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
1426: .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1427: .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1428: -  data - optional location of matrix data.  Set data=NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1429:    to control all matrix memory allocation.

1431:    Output Parameter:
1432: .  A - the matrix

1434:    Notes:
1435:    The dense format is fully compatible with standard Fortran 77
1436:    storage by columns.

1438:    The data input variable is intended primarily for Fortran programmers
1439:    who wish to allocate their own matrix memory space.  Most users should
1440:    set data=NULL (PETSC_NULL_SCALAR for Fortran users).

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

1445:    Level: intermediate

1447: .keywords: matrix,dense, parallel

1449: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1450: @*/
1451: PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
1452: {
1454:   PetscMPIInt    size;

1457:   MatCreate(comm,A);
1458:   MatSetSizes(*A,m,n,M,N);
1459:   MPI_Comm_size(comm,&size);
1460:   if (size > 1) {
1461:     MatSetType(*A,MATMPIDENSE);
1462:     MatMPIDenseSetPreallocation(*A,data);
1463:     if (data) {  /* user provided data array, so no need to assemble */
1464:       MatSetUpMultiply_MPIDense(*A);
1465:       (*A)->assembled = PETSC_TRUE;
1466:     }
1467:   } else {
1468:     MatSetType(*A,MATSEQDENSE);
1469:     MatSeqDenseSetPreallocation(*A,data);
1470:   }
1471:   return(0);
1472: }

1476: static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
1477: {
1478:   Mat            mat;
1479:   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;

1483:   *newmat = 0;
1484:   MatCreate(PetscObjectComm((PetscObject)A),&mat);
1485:   MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
1486:   MatSetType(mat,((PetscObject)A)->type_name);
1487:   a       = (Mat_MPIDense*)mat->data;
1488:   PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));

1490:   mat->factortype   = A->factortype;
1491:   mat->assembled    = PETSC_TRUE;
1492:   mat->preallocated = PETSC_TRUE;

1494:   a->size         = oldmat->size;
1495:   a->rank         = oldmat->rank;
1496:   mat->insertmode = NOT_SET_VALUES;
1497:   a->nvec         = oldmat->nvec;
1498:   a->donotstash   = oldmat->donotstash;

1500:   PetscLayoutReference(A->rmap,&mat->rmap);
1501:   PetscLayoutReference(A->cmap,&mat->cmap);

1503:   MatSetUpMultiply_MPIDense(mat);
1504:   MatDuplicate(oldmat->A,cpvalues,&a->A);
1505:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);

1507:   *newmat = mat;
1508:   return(0);
1509: }

1513: PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat)
1514: {
1516:   PetscMPIInt    rank,size;
1517:   const PetscInt *rowners;
1518:   PetscInt       i,m,n,nz,j,mMax;
1519:   PetscScalar    *array,*vals,*vals_ptr;
1520:   Mat_MPIDense   *a = (Mat_MPIDense*)newmat->data;

1523:   MPI_Comm_rank(comm,&rank);
1524:   MPI_Comm_size(comm,&size);

1526:   /* determine ownership of rows and columns */
1527:   m = (newmat->rmap->n < 0) ? PETSC_DECIDE : newmat->rmap->n;
1528:   n = (newmat->cmap->n < 0) ? PETSC_DECIDE : newmat->cmap->n;

1530:   MatSetSizes(newmat,m,n,M,N);
1531:   if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) {
1532:     MatMPIDenseSetPreallocation(newmat,NULL);
1533:   }
1534:   MatDenseGetArray(newmat,&array);
1535:   MatGetLocalSize(newmat,&m,NULL);
1536:   MatGetOwnershipRanges(newmat,&rowners);
1537:   MPI_Reduce(&m,&mMax,1,MPIU_INT,MPI_MAX,0,comm);
1538:   if (!rank) {
1539:     PetscMalloc1(mMax*N,&vals);

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

1544:     /* insert into matrix-by row (this is why cannot directly read into array */
1545:     vals_ptr = vals;
1546:     for (i=0; i<m; i++) {
1547:       for (j=0; j<N; j++) {
1548:         array[i + j*m] = *vals_ptr++;
1549:       }
1550:     }

1552:     /* read in other processors and ship out */
1553:     for (i=1; i<size; i++) {
1554:       nz   = (rowners[i+1] - rowners[i])*N;
1555:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1556:       MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);
1557:     }
1558:   } else {
1559:     /* receive numeric values */
1560:     PetscMalloc1(m*N,&vals);

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

1565:     /* insert into matrix-by row (this is why cannot directly read into array */
1566:     vals_ptr = vals;
1567:     for (i=0; i<m; i++) {
1568:       for (j=0; j<N; j++) {
1569:         array[i + j*m] = *vals_ptr++;
1570:       }
1571:     }
1572:   }
1573:   MatDenseRestoreArray(newmat,&array);
1574:   PetscFree(vals);
1575:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1576:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1577:   return(0);
1578: }

1582: PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer)
1583: {
1584:   Mat_MPIDense   *a;
1585:   PetscScalar    *vals,*svals;
1586:   MPI_Comm       comm;
1587:   MPI_Status     status;
1588:   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,n,maxnz;
1589:   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
1590:   PetscInt       *ourlens,*procsnz = 0,jj,*mycols,*smycols;
1591:   PetscInt       i,nz,j,rstart,rend;
1592:   int            fd;

1596:   /* force binary viewer to load .info file if it has not yet done so */
1597:   PetscViewerSetUp(viewer);
1598:   PetscObjectGetComm((PetscObject)viewer,&comm);
1599:   MPI_Comm_size(comm,&size);
1600:   MPI_Comm_rank(comm,&rank);
1601:   PetscViewerBinaryGetDescriptor(viewer,&fd);
1602:   if (!rank) {
1603:     PetscBinaryRead(fd,(char*)header,4,PETSC_INT);
1604:     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
1605:   }
1606:   MPI_Bcast(header+1,3,MPIU_INT,0,comm);
1607:   M    = header[1]; N = header[2]; nz = header[3];

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

1613:   if (newmat->rmap->N != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%D) and input matrix has (%D)",M,newmat->rmap->N);
1614:   if (newmat->cmap->N != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%D) and input matrix has (%D)",N,newmat->cmap->N);

1616:   /*
1617:        Handle case where matrix is stored on disk as a dense matrix
1618:   */
1619:   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1620:     MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);
1621:     return(0);
1622:   }

1624:   /* determine ownership of all rows */
1625:   if (newmat->rmap->n < 0) {
1626:     PetscMPIIntCast(M/size + ((M % size) > rank),&m);
1627:   } else {
1628:     PetscMPIIntCast(newmat->rmap->n,&m);
1629:   }
1630:   if (newmat->cmap->n < 0) {
1631:     n = PETSC_DECIDE;
1632:   } else {
1633:     PetscMPIIntCast(newmat->cmap->n,&n);
1634:   }

1636:   PetscMalloc1(size+2,&rowners);
1637:   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1638:   rowners[0] = 0;
1639:   for (i=2; i<=size; i++) {
1640:     rowners[i] += rowners[i-1];
1641:   }
1642:   rstart = rowners[rank];
1643:   rend   = rowners[rank+1];

1645:   /* distribute row lengths to all processors */
1646:   PetscMalloc1(rend-rstart,&ourlens);
1647:   if (!rank) {
1648:     PetscMalloc1(M,&rowlengths);
1649:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
1650:     PetscMalloc1(size,&sndcounts);
1651:     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1652:     MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1653:     PetscFree(sndcounts);
1654:   } else {
1655:     MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);
1656:   }

1658:   if (!rank) {
1659:     /* calculate the number of nonzeros on each processor */
1660:     PetscMalloc1(size,&procsnz);
1661:     PetscMemzero(procsnz,size*sizeof(PetscInt));
1662:     for (i=0; i<size; i++) {
1663:       for (j=rowners[i]; j< rowners[i+1]; j++) {
1664:         procsnz[i] += rowlengths[j];
1665:       }
1666:     }
1667:     PetscFree(rowlengths);

1669:     /* determine max buffer needed and allocate it */
1670:     maxnz = 0;
1671:     for (i=0; i<size; i++) {
1672:       maxnz = PetscMax(maxnz,procsnz[i]);
1673:     }
1674:     PetscMalloc1(maxnz,&cols);

1676:     /* read in my part of the matrix column indices  */
1677:     nz   = procsnz[0];
1678:     PetscMalloc1(nz,&mycols);
1679:     PetscBinaryRead(fd,mycols,nz,PETSC_INT);

1681:     /* read in every one elses and ship off */
1682:     for (i=1; i<size; i++) {
1683:       nz   = procsnz[i];
1684:       PetscBinaryRead(fd,cols,nz,PETSC_INT);
1685:       MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
1686:     }
1687:     PetscFree(cols);
1688:   } else {
1689:     /* determine buffer space needed for message */
1690:     nz = 0;
1691:     for (i=0; i<m; i++) {
1692:       nz += ourlens[i];
1693:     }
1694:     PetscMalloc1(nz+1,&mycols);

1696:     /* receive message of column indices*/
1697:     MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
1698:     MPI_Get_count(&status,MPIU_INT,&maxnz);
1699:     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
1700:   }

1702:   MatSetSizes(newmat,m,n,M,N);
1703:   a = (Mat_MPIDense*)newmat->data;
1704:   if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) {
1705:     MatMPIDenseSetPreallocation(newmat,NULL);
1706:   }

1708:   if (!rank) {
1709:     PetscMalloc1(maxnz,&vals);

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

1715:     /* insert into matrix */
1716:     jj      = rstart;
1717:     smycols = mycols;
1718:     svals   = vals;
1719:     for (i=0; i<m; i++) {
1720:       MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1721:       smycols += ourlens[i];
1722:       svals   += ourlens[i];
1723:       jj++;
1724:     }

1726:     /* read in other processors and ship out */
1727:     for (i=1; i<size; i++) {
1728:       nz   = procsnz[i];
1729:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1730:       MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);
1731:     }
1732:     PetscFree(procsnz);
1733:   } else {
1734:     /* receive numeric values */
1735:     PetscMalloc1(nz+1,&vals);

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

1742:     /* insert into matrix */
1743:     jj      = rstart;
1744:     smycols = mycols;
1745:     svals   = vals;
1746:     for (i=0; i<m; i++) {
1747:       MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);
1748:       smycols += ourlens[i];
1749:       svals   += ourlens[i];
1750:       jj++;
1751:     }
1752:   }
1753:   PetscFree(ourlens);
1754:   PetscFree(vals);
1755:   PetscFree(mycols);
1756:   PetscFree(rowners);

1758:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
1759:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
1760:   return(0);
1761: }

1765: PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
1766: {
1767:   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
1768:   Mat            a,b;
1769:   PetscBool      flg;

1773:   a    = matA->A;
1774:   b    = matB->A;
1775:   MatEqual(a,b,&flg);
1776:   MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1777:   return(0);
1778: }

1782: PetscErrorCode MatDestroy_MatTransMatMult_MPIDense_MPIDense(Mat A)
1783: {
1784:   PetscErrorCode        ierr;
1785:   Mat_MPIDense          *a = (Mat_MPIDense*)A->data;
1786:   Mat_TransMatMultDense *atb = a->atbdense;

1789:   PetscFree3(atb->sendbuf,atb->atbarray,atb->recvcounts);
1790:   (atb->destroy)(A);
1791:   PetscFree(atb);
1792:   return(0);
1793: }

1797: PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C)
1798: {
1799:   Mat_MPIDense   *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data;
1800:   Mat_SeqDense   *aseq=(Mat_SeqDense*)(a->A)->data, *bseq=(Mat_SeqDense*)(b->A)->data;
1801:   Mat_TransMatMultDense *atb = c->atbdense;
1803:   MPI_Comm       comm;
1804:   PetscMPIInt    rank,size,*recvcounts=atb->recvcounts;
1805:   PetscScalar    *carray,*atbarray=atb->atbarray,*sendbuf=atb->sendbuf;
1806:   PetscInt       i,cN=C->cmap->N,cM=C->rmap->N,proc,k,j;
1807:   PetscScalar    _DOne=1.0,_DZero=0.0;
1808:   PetscBLASInt   am,an,bn,aN;
1809:   const PetscInt *ranges;

1812:   PetscObjectGetComm((PetscObject)A,&comm);
1813:   MPI_Comm_rank(comm,&rank);
1814:   MPI_Comm_size(comm,&size);

1816:   /* compute atbarray = aseq^T * bseq */
1817:   PetscBLASIntCast(a->A->cmap->n,&an);
1818:   PetscBLASIntCast(b->A->cmap->n,&bn);
1819:   PetscBLASIntCast(a->A->rmap->n,&am);
1820:   PetscBLASIntCast(A->cmap->N,&aN);
1821:   PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&an,&bn,&am,&_DOne,aseq->v,&aseq->lda,bseq->v,&bseq->lda,&_DZero,atbarray,&aN));
1822: 
1823:   MatGetOwnershipRanges(C,&ranges);
1824:   for (i=0; i<size; i++) recvcounts[i] = (ranges[i+1] - ranges[i])*cN;
1825: 
1826:   /* arrange atbarray into sendbuf */
1827:   k = 0;
1828:   for (proc=0; proc<size; proc++) {
1829:     for (j=0; j<cN; j++) {
1830:       for (i=ranges[proc]; i<ranges[proc+1]; i++) sendbuf[k++] = atbarray[i+j*cM];
1831:     }
1832:   }
1833:   /* sum all atbarray to local values of C */
1834:   MatDenseGetArray(c->A,&carray);
1835:   MPI_Reduce_scatter(sendbuf,carray,recvcounts,MPIU_SCALAR,MPIU_SUM,comm);
1836:   MatDenseRestoreArray(c->A,&carray);
1837:   return(0);
1838: }

1842: PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
1843: {
1844:   PetscErrorCode        ierr;
1845:   Mat                   Cdense;
1846:   MPI_Comm              comm;
1847:   PetscMPIInt           size;
1848:   PetscInt              cm=A->cmap->n,cM,cN=B->cmap->N;
1849:   Mat_MPIDense          *c;
1850:   Mat_TransMatMultDense *atb;

1853:   PetscObjectGetComm((PetscObject)A,&comm);
1854:   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) {
1855:     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
1856:   }

1858:   /* create matrix product Cdense */
1859:   MatCreate(comm,&Cdense);
1860:   MatSetSizes(Cdense,cm,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
1861:   MatSetType(Cdense,MATMPIDENSE);
1862:   MatMPIDenseSetPreallocation(Cdense,NULL);
1863:   MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);
1864:   MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);
1865:   *C   = Cdense;

1867:   /* create data structure for reuse Cdense */
1868:   MPI_Comm_size(comm,&size);
1869:   PetscNew(&atb);
1870:   cM = Cdense->rmap->N;
1871:   PetscMalloc3(cM*cN,&atb->sendbuf,cM*cN,&atb->atbarray,size,&atb->recvcounts);
1872: 
1873:   c                    = (Mat_MPIDense*)Cdense->data;
1874:   c->atbdense          = atb;
1875:   atb->destroy         = Cdense->ops->destroy;
1876:   Cdense->ops->destroy = MatDestroy_MatTransMatMult_MPIDense_MPIDense;
1877:   return(0);
1878: }

1882: PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1883: {

1887:   if (scall == MAT_INITIAL_MATRIX) {
1888:     MatTransposeMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);
1889:   }
1890:   MatTransposeMatMultNumeric_MPIDense_MPIDense(A,B,*C);
1891:   return(0);
1892: }

1896: PetscErrorCode MatDestroy_MatMatMult_MPIDense_MPIDense(Mat A)
1897: {
1898:   PetscErrorCode   ierr;
1899:   Mat_MPIDense     *a = (Mat_MPIDense*)A->data;
1900:   Mat_MatMultDense *ab = a->abdense;

1903:   MatDestroy(&ab->Ce);
1904:   MatDestroy(&ab->Ae);
1905:   MatDestroy(&ab->Be);

1907:   (ab->destroy)(A);
1908:   PetscFree(ab);
1909:   return(0);
1910: }

1912: #if defined(PETSC_HAVE_ELEMENTAL)
1915: PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C)
1916: {
1917:   PetscErrorCode   ierr;
1918:   Mat_MPIDense     *c=(Mat_MPIDense*)C->data;
1919:   Mat_MatMultDense *ab=c->abdense;

1922:   MatConvert_MPIDense_Elemental(A,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Ae);
1923:   MatConvert_MPIDense_Elemental(B,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Be);
1924:   MatMatMultNumeric(ab->Ae,ab->Be,ab->Ce);
1925:   MatConvert(ab->Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C);
1926:   return(0);
1927: }

1931: PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
1932: {
1933:   PetscErrorCode   ierr;
1934:   Mat              Ae,Be,Ce;
1935:   Mat_MPIDense     *c;
1936:   Mat_MatMultDense *ab;

1939:   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) {
1940:     SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
1941:   }

1943:   /* convert A and B to Elemental matrices Ae and Be */
1944:   MatConvert(A,MATELEMENTAL,MAT_INITIAL_MATRIX, &Ae);
1945:   MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX, &Be);

1947:   /* Ce = Ae*Be */
1948:   MatMatMultSymbolic(Ae,Be,fill,&Ce);
1949:   MatMatMultNumeric(Ae,Be,Ce);
1950: 
1951:   /* convert Ce to C */
1952:   MatConvert(Ce,MATMPIDENSE,MAT_INITIAL_MATRIX,C);

1954:   /* create data structure for reuse Cdense */
1955:   PetscNew(&ab);
1956:   c                  = (Mat_MPIDense*)(*C)->data;
1957:   c->abdense         = ab;

1959:   ab->Ae             = Ae;
1960:   ab->Be             = Be;
1961:   ab->Ce             = Ce;
1962:   ab->destroy        = (*C)->ops->destroy;
1963:   (*C)->ops->destroy        = MatDestroy_MatMatMult_MPIDense_MPIDense;
1964:   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIDense_MPIDense;
1965:   return(0);
1966: }

1970: PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1971: {

1975:   if (scall == MAT_INITIAL_MATRIX) { /* simbolic product includes numeric product */
1976:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
1977:     MatMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);
1978:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
1979:   } else {
1980:     PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
1981:     MatMatMultNumeric_MPIDense_MPIDense(A,B,*C);
1982:     PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
1983:   }
1984:   return(0);
1985: }
1986: #endif