Actual source code: mpidense.c

petsc-master 2020-07-09
Report Typos and Errors

  2: /*
  3:    Basic functions for basic parallel dense matrices.
  4: */

  6:  #include <../src/mat/impls/dense/mpi/mpidense.h>
  7:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
  8:  #include <petscblaslapack.h>

 10: /*@

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

 15:     Input Parameter:
 16: .      A - the Seq or MPI dense matrix

 18:     Output Parameter:
 19: .      B - the inner matrix

 21:     Level: intermediate

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

 33:   PetscObjectBaseTypeCompare((PetscObject)A,MATMPIDENSE,&flg);
 34:   if (flg) *B = mat->A;
 35:   else {
 36:     PetscObjectBaseTypeCompare((PetscObject)A,MATSEQDENSE,&flg);
 37:     if (!flg) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)A)->type_name);
 38:     *B = A;
 39:   }
 40:   return(0);
 41: }

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

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

 56: PetscErrorCode MatRestoreRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
 57: {
 58:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
 60:   PetscInt       lrow,rstart = A->rmap->rstart,rend = A->rmap->rend;

 63:   if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows");
 64:   lrow = row - rstart;
 65:   MatRestoreRow(mat->A,lrow,nz,(const PetscInt**)idx,(const PetscScalar**)v);
 66:   return(0);
 67: }

 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:   PetscBool      flg;
 77:   Mat            B;

 80:   MatHasCongruentLayouts(A,&flg);
 81:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only square matrices supported.");
 82:   PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);
 83:   if (!B) { /* This should use MatDenseGetSubMatrix (not create), but we would need a call like MatRestoreDiagonalBlock */

 85:     PetscObjectTypeCompare((PetscObject)mdn->A,MATSEQDENSECUDA,&flg);
 86:     if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for %s. Send an email to petsc-dev@mcs.anl.gov to request this feature",MATSEQDENSECUDA);
 87:     PetscObjectGetComm((PetscObject)(mdn->A),&comm);
 88:     MatCreate(comm,&B);
 89:     MatSetSizes(B,m,m,m,m);
 90:     MatSetType(B,((PetscObject)mdn->A)->type_name);
 91:     MatDenseGetArrayRead(mdn->A,(const PetscScalar**)&array);
 92:     MatSeqDenseSetPreallocation(B,array+m*rstart);
 93:     MatDenseRestoreArrayRead(mdn->A,(const PetscScalar**)&array);
 94:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 95:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 96:     PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);
 97:     *a   = B;
 98:     MatDestroy(&B);
 99:   } else *a = B;
100:   return(0);
101: }

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

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

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: }

159: static PetscErrorCode MatDenseGetLDA_MPIDense(Mat A,PetscInt *lda)
160: {
161:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

165:   MatDenseGetLDA(a->A,lda);
166:   return(0);
167: }

169: static PetscErrorCode MatDenseSetLDA_MPIDense(Mat A,PetscInt lda)
170: {
171:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
172:   PetscBool      iscuda;

176:   if (!a->A) {
177:     if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
178:     PetscLayoutSetUp(A->rmap);
179:     PetscLayoutSetUp(A->cmap);
180:     MatCreate(PETSC_COMM_SELF,&a->A);
181:     PetscLogObjectParent((PetscObject)A,(PetscObject)a->A);
182:     MatSetSizes(a->A,A->rmap->n,A->cmap->N,A->rmap->n,A->cmap->N);
183:     PetscObjectTypeCompare((PetscObject)A,MATMPIDENSECUDA,&iscuda);
184:     MatSetType(a->A,iscuda ? MATSEQDENSECUDA : MATSEQDENSE);
185:   }
186:   MatDenseSetLDA(a->A,lda);
187:   return(0);
188: }

190: static PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar **array)
191: {
192:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

196:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
197:   MatDenseGetArray(a->A,array);
198:   return(0);
199: }

201: static PetscErrorCode MatDenseGetArrayRead_MPIDense(Mat A,const PetscScalar **array)
202: {
203:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

207:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
208:   MatDenseGetArrayRead(a->A,array);
209:   return(0);
210: }

212: static PetscErrorCode MatDenseGetArrayWrite_MPIDense(Mat A,PetscScalar **array)
213: {
214:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

218:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
219:   MatDenseGetArrayWrite(a->A,array);
220:   return(0);
221: }

223: static PetscErrorCode MatDensePlaceArray_MPIDense(Mat A,const PetscScalar *array)
224: {
225:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

229:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
230:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
231:   MatDensePlaceArray(a->A,array);
232:   return(0);
233: }

235: static PetscErrorCode MatDenseResetArray_MPIDense(Mat A)
236: {
237:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

241:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
242:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
243:   MatDenseResetArray(a->A);
244:   return(0);
245: }

247: static PetscErrorCode MatDenseReplaceArray_MPIDense(Mat A,const PetscScalar *array)
248: {
249:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

253:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
254:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
255:   MatDenseReplaceArray(a->A,array);
256:   return(0);
257: }

259: static PetscErrorCode MatCreateSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
260: {
261:   Mat_MPIDense      *mat  = (Mat_MPIDense*)A->data,*newmatd;
262:   PetscErrorCode    ierr;
263:   PetscInt          lda,i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols;
264:   const PetscInt    *irow,*icol;
265:   const PetscScalar *v;
266:   PetscScalar       *bv;
267:   Mat               newmat;
268:   IS                iscol_local;
269:   MPI_Comm          comm_is,comm_mat;

272:   PetscObjectGetComm((PetscObject)A,&comm_mat);
273:   PetscObjectGetComm((PetscObject)iscol,&comm_is);
274:   if (comm_mat != comm_is) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"IS communicator must match matrix communicator");

276:   ISAllGather(iscol,&iscol_local);
277:   ISGetIndices(isrow,&irow);
278:   ISGetIndices(iscol_local,&icol);
279:   ISGetLocalSize(isrow,&nrows);
280:   ISGetLocalSize(iscol,&ncols);
281:   ISGetSize(iscol,&Ncols); /* global number of columns, size of iscol_local */

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

287:   MatGetLocalSize(A,&nlrows,&nlcols);
288:   MatGetOwnershipRange(A,&rstart,&rend);

290:   /* Check submatrix call */
291:   if (scall == MAT_REUSE_MATRIX) {
292:     /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
293:     /* Really need to test rows and column sizes! */
294:     newmat = *B;
295:   } else {
296:     /* Create and fill new matrix */
297:     MatCreate(PetscObjectComm((PetscObject)A),&newmat);
298:     MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);
299:     MatSetType(newmat,((PetscObject)A)->type_name);
300:     MatMPIDenseSetPreallocation(newmat,NULL);
301:   }

303:   /* Now extract the data pointers and do the copy, column at a time */
304:   newmatd = (Mat_MPIDense*)newmat->data;
305:   MatDenseGetArray(newmatd->A,&bv);
306:   MatDenseGetArrayRead(mat->A,&v);
307:   MatDenseGetLDA(mat->A,&lda);
308:   for (i=0; i<Ncols; i++) {
309:     const PetscScalar *av = v + lda*icol[i];
310:     for (j=0; j<nrows; j++) {
311:       *bv++ = av[irow[j] - rstart];
312:     }
313:   }
314:   MatDenseRestoreArrayRead(mat->A,&v);
315:   MatDenseRestoreArray(newmatd->A,&bv);

317:   /* Assemble the matrices so that the correct flags are set */
318:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
319:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);

321:   /* Free work space */
322:   ISRestoreIndices(isrow,&irow);
323:   ISRestoreIndices(iscol_local,&icol);
324:   ISDestroy(&iscol_local);
325:   *B   = newmat;
326:   return(0);
327: }

329: PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar **array)
330: {
331:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

335:   MatDenseRestoreArray(a->A,array);
336:   return(0);
337: }

339: PetscErrorCode MatDenseRestoreArrayRead_MPIDense(Mat A,const PetscScalar **array)
340: {
341:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

345:   MatDenseRestoreArrayRead(a->A,array);
346:   return(0);
347: }

349: PetscErrorCode MatDenseRestoreArrayWrite_MPIDense(Mat A,PetscScalar **array)
350: {
351:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

355:   MatDenseRestoreArrayWrite(a->A,array);
356:   return(0);
357: }

359: PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
360: {
361:   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
363:   PetscInt       nstash,reallocs;

366:   if (mdn->donotstash || mat->nooffprocentries) return(0);

368:   MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
369:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
370:   PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
371:   return(0);
372: }

374: PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
375: {
376:   Mat_MPIDense   *mdn=(Mat_MPIDense*)mat->data;
378:   PetscInt       i,*row,*col,flg,j,rstart,ncols;
379:   PetscMPIInt    n;
380:   PetscScalar    *val;

383:   if (!mdn->donotstash && !mat->nooffprocentries) {
384:     /*  wait on receives */
385:     while (1) {
386:       MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
387:       if (!flg) break;

389:       for (i=0; i<n;) {
390:         /* Now identify the consecutive vals belonging to the same row */
391:         for (j=i,rstart=row[j]; j<n; j++) {
392:           if (row[j] != rstart) break;
393:         }
394:         if (j < n) ncols = j-i;
395:         else       ncols = n-i;
396:         /* Now assemble all these values with a single function call */
397:         MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);
398:         i    = j;
399:       }
400:     }
401:     MatStashScatterEnd_Private(&mat->stash);
402:   }

404:   MatAssemblyBegin(mdn->A,mode);
405:   MatAssemblyEnd(mdn->A,mode);

407:   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
408:     MatSetUpMultiply_MPIDense(mat);
409:   }
410:   return(0);
411: }

413: PetscErrorCode MatZeroEntries_MPIDense(Mat A)
414: {
416:   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;

419:   MatZeroEntries(l->A);
420:   return(0);
421: }

423: PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
424: {
425:   Mat_MPIDense      *l = (Mat_MPIDense*)A->data;
426:   PetscErrorCode    ierr;
427:   PetscInt          i,len,*lrows;

430:   /* get locally owned rows */
431:   PetscLayoutMapLocal(A->rmap,n,rows,&len,&lrows,NULL);
432:   /* fix right hand side if needed */
433:   if (x && b) {
434:     const PetscScalar *xx;
435:     PetscScalar       *bb;

437:     VecGetArrayRead(x, &xx);
438:     VecGetArrayWrite(b, &bb);
439:     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
440:     VecRestoreArrayRead(x, &xx);
441:     VecRestoreArrayWrite(b, &bb);
442:   }
443:   MatZeroRows(l->A,len,lrows,0.0,NULL,NULL);
444:   if (diag != 0.0) {
445:     Vec d;

447:     MatCreateVecs(A,NULL,&d);
448:     VecSet(d,diag);
449:     MatDiagonalSet(A,d,INSERT_VALUES);
450:     VecDestroy(&d);
451:   }
452:   PetscFree(lrows);
453:   return(0);
454: }

456: PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
457: PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
458: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
459: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);

461: PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
462: {
463:   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
464:   PetscErrorCode    ierr;
465:   const PetscScalar *ax;
466:   PetscScalar       *ay;

469:   VecGetArrayReadInPlace(xx,&ax);
470:   VecGetArrayInPlace(mdn->lvec,&ay);
471:   PetscSFBcastBegin(mdn->Mvctx,MPIU_SCALAR,ax,ay);
472:   PetscSFBcastEnd(mdn->Mvctx,MPIU_SCALAR,ax,ay);
473:   VecRestoreArrayInPlace(mdn->lvec,&ay);
474:   VecRestoreArrayReadInPlace(xx,&ax);
475:   (*mdn->A->ops->mult)(mdn->A,mdn->lvec,yy);
476:   return(0);
477: }

479: PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
480: {
481:   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
482:   PetscErrorCode    ierr;
483:   const PetscScalar *ax;
484:   PetscScalar       *ay;

487:   VecGetArrayReadInPlace(xx,&ax);
488:   VecGetArrayInPlace(mdn->lvec,&ay);
489:   PetscSFBcastBegin(mdn->Mvctx,MPIU_SCALAR,ax,ay);
490:   PetscSFBcastEnd(mdn->Mvctx,MPIU_SCALAR,ax,ay);
491:   VecRestoreArrayInPlace(mdn->lvec,&ay);
492:   VecRestoreArrayReadInPlace(xx,&ax);
493:   (*mdn->A->ops->multadd)(mdn->A,mdn->lvec,yy,zz);
494:   return(0);
495: }

497: PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
498: {
499:   Mat_MPIDense      *a = (Mat_MPIDense*)A->data;
500:   PetscErrorCode    ierr;
501:   const PetscScalar *ax;
502:   PetscScalar       *ay;

505:   VecSet(yy,0.0);
506:   (*a->A->ops->multtranspose)(a->A,xx,a->lvec);
507:   VecGetArrayReadInPlace(a->lvec,&ax);
508:   VecGetArrayInPlace(yy,&ay);
509:   PetscSFReduceBegin(a->Mvctx,MPIU_SCALAR,ax,ay,MPIU_SUM);
510:   PetscSFReduceEnd(a->Mvctx,MPIU_SCALAR,ax,ay,MPIU_SUM);
511:   VecRestoreArrayReadInPlace(a->lvec,&ax);
512:   VecRestoreArrayInPlace(yy,&ay);
513:   return(0);
514: }

516: PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
517: {
518:   Mat_MPIDense      *a = (Mat_MPIDense*)A->data;
519:   PetscErrorCode    ierr;
520:   const PetscScalar *ax;
521:   PetscScalar       *ay;

524:   VecCopy(yy,zz);
525:   (*a->A->ops->multtranspose)(a->A,xx,a->lvec);
526:   VecGetArrayReadInPlace(a->lvec,&ax);
527:   VecGetArrayInPlace(zz,&ay);
528:   PetscSFReduceBegin(a->Mvctx,MPIU_SCALAR,ax,ay,MPIU_SUM);
529:   PetscSFReduceEnd(a->Mvctx,MPIU_SCALAR,ax,ay,MPIU_SUM);
530:   VecRestoreArrayReadInPlace(a->lvec,&ax);
531:   VecRestoreArrayInPlace(zz,&ay);
532:   return(0);
533: }

535: PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
536: {
537:   Mat_MPIDense      *a    = (Mat_MPIDense*)A->data;
538:   PetscErrorCode    ierr;
539:   PetscInt          lda,len,i,n,m = A->rmap->n,radd;
540:   PetscScalar       *x,zero = 0.0;
541:   const PetscScalar *av;

544:   VecSet(v,zero);
545:   VecGetArray(v,&x);
546:   VecGetSize(v,&n);
547:   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
548:   len  = PetscMin(a->A->rmap->n,a->A->cmap->n);
549:   radd = A->rmap->rstart*m;
550:   MatDenseGetArrayRead(a->A,&av);
551:   MatDenseGetLDA(a->A,&lda);
552:   for (i=0; i<len; i++) {
553:     x[i] = av[radd + i*lda + i];
554:   }
555:   MatDenseRestoreArrayRead(a->A,&av);
556:   VecRestoreArray(v,&x);
557:   return(0);
558: }

560: PetscErrorCode MatDestroy_MPIDense(Mat mat)
561: {
562:   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;

566: #if defined(PETSC_USE_LOG)
567:   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N);
568: #endif
569:   MatStashDestroy_Private(&mat->stash);
570:   if (mdn->vecinuse) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
571:   if (mdn->matinuse) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
572:   MatDestroy(&mdn->A);
573:   VecDestroy(&mdn->lvec);
574:   PetscSFDestroy(&mdn->Mvctx);
575:   VecDestroy(&mdn->cvec);
576:   MatDestroy(&mdn->cmat);

578:   PetscFree(mat->data);
579:   PetscObjectChangeTypeName((PetscObject)mat,0);

581:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",NULL);
582:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseSetLDA_C",NULL);
583:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);
584:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);
585:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",NULL);
586:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",NULL);
587:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayWrite_C",NULL);
588:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayWrite_C",NULL);
589:   PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);
590:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);
591:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseReplaceArray_C",NULL);
592: #if defined(PETSC_HAVE_ELEMENTAL)
593:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",NULL);
594: #endif
595: #if defined(PETSC_HAVE_SCALAPACK)
596:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_scalapack_C",NULL);
597: #endif
598:   PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",NULL);
599:   PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_mpiaij_mpidense_C",NULL);
600:   PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_mpidense_mpiaij_C",NULL);
601: #if defined (PETSC_HAVE_CUDA)
602:   PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_mpiaijcusparse_mpidense_C",NULL);
603:   PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_mpidense_mpiaijcusparse_C",NULL);
604: #endif
605:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL);
606:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL);
607:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVec_C",NULL);
608:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVec_C",NULL);
609:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecRead_C",NULL);
610:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecRead_C",NULL);
611:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecWrite_C",NULL);
612:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecWrite_C",NULL);
613:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetSubMatrix_C",NULL);
614:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreSubMatrix_C",NULL);
615:   return(0);
616: }

618: PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer);

620:  #include <petscdraw.h>
621: static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
622: {
623:   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
624:   PetscErrorCode    ierr;
625:   PetscMPIInt       rank;
626:   PetscViewerType   vtype;
627:   PetscBool         iascii,isdraw;
628:   PetscViewer       sviewer;
629:   PetscViewerFormat format;

632:   MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
633:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
634:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
635:   if (iascii) {
636:     PetscViewerGetType(viewer,&vtype);
637:     PetscViewerGetFormat(viewer,&format);
638:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
639:       MatInfo info;
640:       MatGetInfo(mat,MAT_LOCAL,&info);
641:       PetscViewerASCIIPushSynchronized(viewer);
642:       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);
643:       PetscViewerFlush(viewer);
644:       PetscViewerASCIIPopSynchronized(viewer);
645:       PetscSFView(mdn->Mvctx,viewer);
646:       return(0);
647:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
648:       return(0);
649:     }
650:   } else if (isdraw) {
651:     PetscDraw draw;
652:     PetscBool isnull;

654:     PetscViewerDrawGetDraw(viewer,0,&draw);
655:     PetscDrawIsNull(draw,&isnull);
656:     if (isnull) return(0);
657:   }

659:   {
660:     /* assemble the entire matrix onto first processor. */
661:     Mat         A;
662:     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
663:     PetscInt    *cols;
664:     PetscScalar *vals;

666:     MatCreate(PetscObjectComm((PetscObject)mat),&A);
667:     if (!rank) {
668:       MatSetSizes(A,M,N,M,N);
669:     } else {
670:       MatSetSizes(A,0,0,M,N);
671:     }
672:     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
673:     MatSetType(A,MATMPIDENSE);
674:     MatMPIDenseSetPreallocation(A,NULL);
675:     PetscLogObjectParent((PetscObject)mat,(PetscObject)A);

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

681:     row = mat->rmap->rstart;
682:     m   = mdn->A->rmap->n;
683:     for (i=0; i<m; i++) {
684:       MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);
685:       MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);
686:       MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);
687:       row++;
688:     }

690:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
691:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
692:     PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
693:     if (!rank) {
694:       PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);
695:       MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);
696:     }
697:     PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
698:     PetscViewerFlush(viewer);
699:     MatDestroy(&A);
700:   }
701:   return(0);
702: }

704: PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
705: {
707:   PetscBool      iascii,isbinary,isdraw,issocket;

710:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
711:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
712:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
713:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);

715:   if (iascii || issocket || isdraw) {
716:     MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);
717:   } else if (isbinary) {
718:     MatView_Dense_Binary(mat,viewer);
719:   }
720:   return(0);
721: }

723: PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
724: {
725:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
726:   Mat            mdn  = mat->A;
728:   PetscLogDouble isend[5],irecv[5];

731:   info->block_size = 1.0;

733:   MatGetInfo(mdn,MAT_LOCAL,info);

735:   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
736:   isend[3] = info->memory;  isend[4] = info->mallocs;
737:   if (flag == MAT_LOCAL) {
738:     info->nz_used      = isend[0];
739:     info->nz_allocated = isend[1];
740:     info->nz_unneeded  = isend[2];
741:     info->memory       = isend[3];
742:     info->mallocs      = isend[4];
743:   } else if (flag == MAT_GLOBAL_MAX) {
744:     MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_MAX,PetscObjectComm((PetscObject)A));

746:     info->nz_used      = irecv[0];
747:     info->nz_allocated = irecv[1];
748:     info->nz_unneeded  = irecv[2];
749:     info->memory       = irecv[3];
750:     info->mallocs      = irecv[4];
751:   } else if (flag == MAT_GLOBAL_SUM) {
752:     MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_SUM,PetscObjectComm((PetscObject)A));

754:     info->nz_used      = irecv[0];
755:     info->nz_allocated = irecv[1];
756:     info->nz_unneeded  = irecv[2];
757:     info->memory       = irecv[3];
758:     info->mallocs      = irecv[4];
759:   }
760:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
761:   info->fill_ratio_needed = 0;
762:   info->factor_mallocs    = 0;
763:   return(0);
764: }

766: PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg)
767: {
768:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

772:   switch (op) {
773:   case MAT_NEW_NONZERO_LOCATIONS:
774:   case MAT_NEW_NONZERO_LOCATION_ERR:
775:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
776:     MatCheckPreallocated(A,1);
777:     MatSetOption(a->A,op,flg);
778:     break;
779:   case MAT_ROW_ORIENTED:
780:     MatCheckPreallocated(A,1);
781:     a->roworiented = flg;
782:     MatSetOption(a->A,op,flg);
783:     break;
784:   case MAT_NEW_DIAGONALS:
785:   case MAT_KEEP_NONZERO_PATTERN:
786:   case MAT_USE_HASH_TABLE:
787:   case MAT_SORTED_FULL:
788:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
789:     break;
790:   case MAT_IGNORE_OFF_PROC_ENTRIES:
791:     a->donotstash = flg;
792:     break;
793:   case MAT_SYMMETRIC:
794:   case MAT_STRUCTURALLY_SYMMETRIC:
795:   case MAT_HERMITIAN:
796:   case MAT_SYMMETRY_ETERNAL:
797:   case MAT_IGNORE_LOWER_TRIANGULAR:
798:   case MAT_IGNORE_ZERO_ENTRIES:
799:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
800:     break;
801:   default:
802:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
803:   }
804:   return(0);
805: }

807: PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
808: {
809:   Mat_MPIDense      *mdn = (Mat_MPIDense*)A->data;
810:   const PetscScalar *l;
811:   PetscScalar       x,*v,*vv,*r;
812:   PetscErrorCode    ierr;
813:   PetscInt          i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n,lda;

816:   MatDenseGetArray(mdn->A,&vv);
817:   MatDenseGetLDA(mdn->A,&lda);
818:   MatGetLocalSize(A,&s2,&s3);
819:   if (ll) {
820:     VecGetLocalSize(ll,&s2a);
821:     if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %D != %D", s2a, s2);
822:     VecGetArrayRead(ll,&l);
823:     for (i=0; i<m; i++) {
824:       x = l[i];
825:       v = vv + i;
826:       for (j=0; j<n; j++) { (*v) *= x; v+= lda;}
827:     }
828:     VecRestoreArrayRead(ll,&l);
829:     PetscLogFlops(1.0*n*m);
830:   }
831:   if (rr) {
832:     const PetscScalar *ar;

834:     VecGetLocalSize(rr,&s3a);
835:     if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
836:     VecGetArrayRead(rr,&ar);
837:     VecGetArray(mdn->lvec,&r);
838:     PetscSFBcastBegin(mdn->Mvctx,MPIU_SCALAR,ar,r);
839:     PetscSFBcastEnd(mdn->Mvctx,MPIU_SCALAR,ar,r);
840:     VecRestoreArrayRead(rr,&ar);
841:     for (i=0; i<n; i++) {
842:       x = r[i];
843:       v = vv + i*lda;
844:       for (j=0; j<m; j++) (*v++) *= x;
845:     }
846:     VecRestoreArray(mdn->lvec,&r);
847:     PetscLogFlops(1.0*n*m);
848:   }
849:   MatDenseRestoreArray(mdn->A,&vv);
850:   return(0);
851: }

853: PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
854: {
855:   Mat_MPIDense      *mdn = (Mat_MPIDense*)A->data;
856:   PetscErrorCode    ierr;
857:   PetscInt          i,j;
858:   PetscMPIInt       size;
859:   PetscReal         sum = 0.0;
860:   const PetscScalar *av,*v;

863:   MatDenseGetArrayRead(mdn->A,&av);
864:   v    = av;
865:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
866:   if (size == 1) {
867:     MatNorm(mdn->A,type,nrm);
868:   } else {
869:     if (type == NORM_FROBENIUS) {
870:       for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) {
871:         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
872:       }
873:       MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));
874:       *nrm = PetscSqrtReal(*nrm);
875:       PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);
876:     } else if (type == NORM_1) {
877:       PetscReal *tmp,*tmp2;
878:       PetscCalloc2(A->cmap->N,&tmp,A->cmap->N,&tmp2);
879:       *nrm = 0.0;
880:       v    = av;
881:       for (j=0; j<mdn->A->cmap->n; j++) {
882:         for (i=0; i<mdn->A->rmap->n; i++) {
883:           tmp[j] += PetscAbsScalar(*v);  v++;
884:         }
885:       }
886:       MPIU_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));
887:       for (j=0; j<A->cmap->N; j++) {
888:         if (tmp2[j] > *nrm) *nrm = tmp2[j];
889:       }
890:       PetscFree2(tmp,tmp2);
891:       PetscLogFlops(A->cmap->n*A->rmap->n);
892:     } else if (type == NORM_INFINITY) { /* max row norm */
893:       PetscReal ntemp;
894:       MatNorm(mdn->A,type,&ntemp);
895:       MPIU_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));
896:     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm");
897:   }
898:   MatDenseRestoreArrayRead(mdn->A,&av);
899:   return(0);
900: }

902: PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
903: {
904:   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
905:   Mat            B;
906:   PetscInt       M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart;
908:   PetscInt       j,i,lda;
909:   PetscScalar    *v;

912:   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
913:     MatCreate(PetscObjectComm((PetscObject)A),&B);
914:     MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);
915:     MatSetType(B,((PetscObject)A)->type_name);
916:     MatMPIDenseSetPreallocation(B,NULL);
917:   } else B = *matout;

919:   m    = a->A->rmap->n; n = a->A->cmap->n;
920:   MatDenseGetArrayRead(a->A,(const PetscScalar**)&v);
921:   MatDenseGetLDA(a->A,&lda);
922:   PetscMalloc1(m,&rwork);
923:   for (i=0; i<m; i++) rwork[i] = rstart + i;
924:   for (j=0; j<n; j++) {
925:     MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);
926:     v   += lda;
927:   }
928:   MatDenseRestoreArrayRead(a->A,(const PetscScalar**)&v);
929:   PetscFree(rwork);
930:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
931:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
932:   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
933:     *matout = B;
934:   } else {
935:     MatHeaderMerge(A,&B);
936:   }
937:   return(0);
938: }

940: static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*);
941: PETSC_INTERN PetscErrorCode MatScale_MPIDense(Mat,PetscScalar);

943: PetscErrorCode MatSetUp_MPIDense(Mat A)
944: {

948:   PetscLayoutSetUp(A->rmap);
949:   PetscLayoutSetUp(A->cmap);
950:   if (!A->preallocated) {
951:     MatMPIDenseSetPreallocation(A,0);
952:   }
953:   return(0);
954: }

956: PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
957: {
959:   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;

962:   MatAXPY(A->A,alpha,B->A,str);
963:   return(0);
964: }

966: PetscErrorCode MatConjugate_MPIDense(Mat mat)
967: {
968:   Mat_MPIDense   *a = (Mat_MPIDense*)mat->data;

972:   MatConjugate(a->A);
973:   return(0);
974: }

976: PetscErrorCode MatRealPart_MPIDense(Mat A)
977: {
978:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

982:   MatRealPart(a->A);
983:   return(0);
984: }

986: PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
987: {
988:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

992:   MatImaginaryPart(a->A);
993:   return(0);
994: }

996: static PetscErrorCode MatGetColumnVector_MPIDense(Mat A,Vec v,PetscInt col)
997: {
999:   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;

1002:   if (!a->A) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Missing local matrix");
1003:   if (!a->A->ops->getcolumnvector) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Missing get column operation");
1004:   (*a->A->ops->getcolumnvector)(a->A,v,col);
1005:   return(0);
1006: }

1008: PETSC_INTERN PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);

1010: PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
1011: {
1013:   PetscInt       i,n;
1014:   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
1015:   PetscReal      *work;

1018:   MatGetSize(A,NULL,&n);
1019:   PetscMalloc1(n,&work);
1020:   MatGetColumnNorms_SeqDense(a->A,type,work);
1021:   if (type == NORM_2) {
1022:     for (i=0; i<n; i++) work[i] *= work[i];
1023:   }
1024:   if (type == NORM_INFINITY) {
1025:     MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);
1026:   } else {
1027:     MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);
1028:   }
1029:   PetscFree(work);
1030:   if (type == NORM_2) {
1031:     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
1032:   }
1033:   return(0);
1034: }

1036: #if defined(PETSC_HAVE_CUDA)
1037: static PetscErrorCode MatDenseGetColumnVec_MPIDenseCUDA(Mat A,PetscInt col,Vec *v)
1038: {
1039:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1041:   PetscInt       lda;

1044:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1045:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1046:   if (!a->cvec) {
1047:     VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,A->rmap->N,NULL,&a->cvec);
1048:     PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec);
1049:   }
1050:   a->vecinuse = col + 1;
1051:   MatDenseGetLDA(a->A,&lda);
1052:   MatDenseCUDAGetArray(a->A,(PetscScalar**)&a->ptrinuse);
1053:   VecCUDAPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)lda);
1054:   *v   = a->cvec;
1055:   return(0);
1056: }

1058: static PetscErrorCode MatDenseRestoreColumnVec_MPIDenseCUDA(Mat A,PetscInt col,Vec *v)
1059: {
1060:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1064:   if (!a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
1065:   if (!a->cvec) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column vector");
1066:   a->vecinuse = 0;
1067:   MatDenseCUDARestoreArray(a->A,(PetscScalar**)&a->ptrinuse);
1068:   VecCUDAResetArray(a->cvec);
1069:   *v   = NULL;
1070:   return(0);
1071: }

1073: static PetscErrorCode MatDenseGetColumnVecRead_MPIDenseCUDA(Mat A,PetscInt col,Vec *v)
1074: {
1075:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1076:   PetscInt       lda;

1080:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1081:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1082:   if (!a->cvec) {
1083:     VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,A->rmap->N,NULL,&a->cvec);
1084:     PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec);
1085:   }
1086:   a->vecinuse = col + 1;
1087:   MatDenseGetLDA(a->A,&lda);
1088:   MatDenseCUDAGetArrayRead(a->A,&a->ptrinuse);
1089:   VecCUDAPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)lda);
1090:   VecLockReadPush(a->cvec);
1091:   *v   = a->cvec;
1092:   return(0);
1093: }

1095: static PetscErrorCode MatDenseRestoreColumnVecRead_MPIDenseCUDA(Mat A,PetscInt col,Vec *v)
1096: {
1097:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1101:   if (!a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
1102:   if (!a->cvec) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing internal column vector");
1103:   a->vecinuse = 0;
1104:   MatDenseCUDARestoreArrayRead(a->A,&a->ptrinuse);
1105:   VecLockReadPop(a->cvec);
1106:   VecCUDAResetArray(a->cvec);
1107:   *v   = NULL;
1108:   return(0);
1109: }

1111: static PetscErrorCode MatDenseGetColumnVecWrite_MPIDenseCUDA(Mat A,PetscInt col,Vec *v)
1112: {
1113:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1115:   PetscInt       lda;

1118:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1119:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1120:   if (!a->cvec) {
1121:     VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,A->rmap->N,NULL,&a->cvec);
1122:     PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec);
1123:   }
1124:   a->vecinuse = col + 1;
1125:   MatDenseGetLDA(a->A,&lda);
1126:   MatDenseCUDAGetArrayWrite(a->A,(PetscScalar**)&a->ptrinuse);
1127:   VecCUDAPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)lda);
1128:   *v   = a->cvec;
1129:   return(0);
1130: }

1132: static PetscErrorCode MatDenseRestoreColumnVecWrite_MPIDenseCUDA(Mat A,PetscInt col,Vec *v)
1133: {
1134:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1138:   if (!a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
1139:   if (!a->cvec) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing internal column vector");
1140:   a->vecinuse = 0;
1141:   MatDenseCUDARestoreArrayWrite(a->A,(PetscScalar**)&a->ptrinuse);
1142:   VecCUDAResetArray(a->cvec);
1143:   *v   = NULL;
1144:   return(0);
1145: }

1147: static PetscErrorCode MatDenseCUDAPlaceArray_MPIDenseCUDA(Mat A, const PetscScalar *a)
1148: {
1149:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1153:   if (l->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1154:   if (l->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1155:   MatDenseCUDAPlaceArray(l->A,a);
1156:   return(0);
1157: }

1159: static PetscErrorCode MatDenseCUDAResetArray_MPIDenseCUDA(Mat A)
1160: {
1161:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1165:   if (l->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1166:   if (l->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1167:   MatDenseCUDAResetArray(l->A);
1168:   return(0);
1169: }

1171: static PetscErrorCode MatDenseCUDAReplaceArray_MPIDenseCUDA(Mat A, const PetscScalar *a)
1172: {
1173:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1177:   if (l->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1178:   if (l->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1179:   MatDenseCUDAReplaceArray(l->A,a);
1180:   return(0);
1181: }

1183: static PetscErrorCode MatDenseCUDAGetArrayWrite_MPIDenseCUDA(Mat A, PetscScalar **a)
1184: {
1185:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1189:   MatDenseCUDAGetArrayWrite(l->A,a);
1190:   return(0);
1191: }

1193: static PetscErrorCode MatDenseCUDARestoreArrayWrite_MPIDenseCUDA(Mat A, PetscScalar **a)
1194: {
1195:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1199:   MatDenseCUDARestoreArrayWrite(l->A,a);
1200:   return(0);
1201: }

1203: static PetscErrorCode MatDenseCUDAGetArrayRead_MPIDenseCUDA(Mat A, const PetscScalar **a)
1204: {
1205:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1209:   MatDenseCUDAGetArrayRead(l->A,a);
1210:   return(0);
1211: }

1213: static PetscErrorCode MatDenseCUDARestoreArrayRead_MPIDenseCUDA(Mat A, const PetscScalar **a)
1214: {
1215:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1219:   MatDenseCUDARestoreArrayRead(l->A,a);
1220:   return(0);
1221: }

1223: static PetscErrorCode MatDenseCUDAGetArray_MPIDenseCUDA(Mat A, PetscScalar **a)
1224: {
1225:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1229:   MatDenseCUDAGetArray(l->A,a);
1230:   return(0);
1231: }

1233: static PetscErrorCode MatDenseCUDARestoreArray_MPIDenseCUDA(Mat A, PetscScalar **a)
1234: {
1235:   Mat_MPIDense   *l = (Mat_MPIDense*) A->data;

1239:   MatDenseCUDARestoreArray(l->A,a);
1240:   return(0);
1241: }

1243: static PetscErrorCode MatDenseGetColumnVecWrite_MPIDense(Mat,PetscInt,Vec*);
1244: static PetscErrorCode MatDenseGetColumnVecRead_MPIDense(Mat,PetscInt,Vec*);
1245: static PetscErrorCode MatDenseGetColumnVec_MPIDense(Mat,PetscInt,Vec*);
1246: static PetscErrorCode MatDenseRestoreColumnVecWrite_MPIDense(Mat,PetscInt,Vec*);
1247: static PetscErrorCode MatDenseRestoreColumnVecRead_MPIDense(Mat,PetscInt,Vec*);
1248: static PetscErrorCode MatDenseRestoreColumnVec_MPIDense(Mat,PetscInt,Vec*);
1249: static PetscErrorCode MatDenseRestoreSubMatrix_MPIDense(Mat,Mat*);

1251: static PetscErrorCode MatBindToCPU_MPIDenseCUDA(Mat mat,PetscBool bind)
1252: {
1253:   Mat_MPIDense   *d = (Mat_MPIDense*)mat->data;

1257:   if (d->vecinuse) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1258:   if (d->matinuse) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1259:   if (d->A) {
1260:     MatBindToCPU(d->A,bind);
1261:   }
1262:   mat->boundtocpu = bind;
1263:   if (!bind) {
1264:     PetscBool iscuda;

1266:     PetscObjectTypeCompare((PetscObject)d->cvec,VECMPICUDA,&iscuda);
1267:     if (!iscuda) {
1268:       VecDestroy(&d->cvec);
1269:     }
1270:     PetscObjectTypeCompare((PetscObject)d->cmat,MATMPIDENSECUDA,&iscuda);
1271:     if (!iscuda) {
1272:       MatDestroy(&d->cmat);
1273:     }
1274:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVec_C",MatDenseGetColumnVec_MPIDenseCUDA);
1275:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVec_C",MatDenseRestoreColumnVec_MPIDenseCUDA);
1276:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecRead_C",MatDenseGetColumnVecRead_MPIDenseCUDA);
1277:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecRead_C",MatDenseRestoreColumnVecRead_MPIDenseCUDA);
1278:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecWrite_C",MatDenseGetColumnVecWrite_MPIDenseCUDA);
1279:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecWrite_C",MatDenseRestoreColumnVecWrite_MPIDenseCUDA);
1280:   } else {
1281:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVec_C",MatDenseGetColumnVec_MPIDense);
1282:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVec_C",MatDenseRestoreColumnVec_MPIDense);
1283:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecRead_C",MatDenseGetColumnVecRead_MPIDense);
1284:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecRead_C",MatDenseRestoreColumnVecRead_MPIDense);
1285:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecWrite_C",MatDenseGetColumnVecWrite_MPIDense);
1286:     PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecWrite_C",MatDenseRestoreColumnVecWrite_MPIDense);
1287:   }
1288:   if (d->cmat) {
1289:     MatBindToCPU(d->cmat,bind);
1290:   }
1291:   return(0);
1292: }

1294: PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat A, PetscScalar *d_data)
1295: {
1296:   Mat_MPIDense   *d = (Mat_MPIDense*)A->data;
1298:   PetscBool      iscuda;

1302:   PetscObjectTypeCompare((PetscObject)A,MATMPIDENSECUDA,&iscuda);
1303:   if (!iscuda) return(0);
1304:   PetscLayoutSetUp(A->rmap);
1305:   PetscLayoutSetUp(A->cmap);
1306:   if (!d->A) {
1307:     MatCreate(PETSC_COMM_SELF,&d->A);
1308:     PetscLogObjectParent((PetscObject)A,(PetscObject)d->A);
1309:     MatSetSizes(d->A,A->rmap->n,A->cmap->N,A->rmap->n,A->cmap->N);
1310:   }
1311:   MatSetType(d->A,MATSEQDENSECUDA);
1312:   MatSeqDenseCUDASetPreallocation(d->A,d_data);
1313:   A->preallocated = PETSC_TRUE;
1314:   return(0);
1315: }
1316: #endif

1318: static PetscErrorCode MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
1319: {
1320:   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;

1324:   MatSetRandom(d->A,rctx);
1325:   return(0);
1326: }

1328: static PetscErrorCode MatMissingDiagonal_MPIDense(Mat A,PetscBool  *missing,PetscInt *d)
1329: {
1331:   *missing = PETSC_FALSE;
1332:   return(0);
1333: }

1335: static PetscErrorCode MatMatTransposeMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat);
1336: static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
1337: static PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
1338: static PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat);
1339: static PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*);
1340: static PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer);

1342: /* -------------------------------------------------------------------*/
1343: static struct _MatOps MatOps_Values = { MatSetValues_MPIDense,
1344:                                         MatGetRow_MPIDense,
1345:                                         MatRestoreRow_MPIDense,
1346:                                         MatMult_MPIDense,
1347:                                 /*  4*/ MatMultAdd_MPIDense,
1348:                                         MatMultTranspose_MPIDense,
1349:                                         MatMultTransposeAdd_MPIDense,
1350:                                         0,
1351:                                         0,
1352:                                         0,
1353:                                 /* 10*/ 0,
1354:                                         0,
1355:                                         0,
1356:                                         0,
1357:                                         MatTranspose_MPIDense,
1358:                                 /* 15*/ MatGetInfo_MPIDense,
1359:                                         MatEqual_MPIDense,
1360:                                         MatGetDiagonal_MPIDense,
1361:                                         MatDiagonalScale_MPIDense,
1362:                                         MatNorm_MPIDense,
1363:                                 /* 20*/ MatAssemblyBegin_MPIDense,
1364:                                         MatAssemblyEnd_MPIDense,
1365:                                         MatSetOption_MPIDense,
1366:                                         MatZeroEntries_MPIDense,
1367:                                 /* 24*/ MatZeroRows_MPIDense,
1368:                                         0,
1369:                                         0,
1370:                                         0,
1371:                                         0,
1372:                                 /* 29*/ MatSetUp_MPIDense,
1373:                                         0,
1374:                                         0,
1375:                                         MatGetDiagonalBlock_MPIDense,
1376:                                         0,
1377:                                 /* 34*/ MatDuplicate_MPIDense,
1378:                                         0,
1379:                                         0,
1380:                                         0,
1381:                                         0,
1382:                                 /* 39*/ MatAXPY_MPIDense,
1383:                                         MatCreateSubMatrices_MPIDense,
1384:                                         0,
1385:                                         MatGetValues_MPIDense,
1386:                                         0,
1387:                                 /* 44*/ 0,
1388:                                         MatScale_MPIDense,
1389:                                         MatShift_Basic,
1390:                                         0,
1391:                                         0,
1392:                                 /* 49*/ MatSetRandom_MPIDense,
1393:                                         0,
1394:                                         0,
1395:                                         0,
1396:                                         0,
1397:                                 /* 54*/ 0,
1398:                                         0,
1399:                                         0,
1400:                                         0,
1401:                                         0,
1402:                                 /* 59*/ MatCreateSubMatrix_MPIDense,
1403:                                         MatDestroy_MPIDense,
1404:                                         MatView_MPIDense,
1405:                                         0,
1406:                                         0,
1407:                                 /* 64*/ 0,
1408:                                         0,
1409:                                         0,
1410:                                         0,
1411:                                         0,
1412:                                 /* 69*/ 0,
1413:                                         0,
1414:                                         0,
1415:                                         0,
1416:                                         0,
1417:                                 /* 74*/ 0,
1418:                                         0,
1419:                                         0,
1420:                                         0,
1421:                                         0,
1422:                                 /* 79*/ 0,
1423:                                         0,
1424:                                         0,
1425:                                         0,
1426:                                 /* 83*/ MatLoad_MPIDense,
1427:                                         0,
1428:                                         0,
1429:                                         0,
1430:                                         0,
1431:                                         0,
1432:                                 /* 89*/ 0,
1433:                                         0,
1434:                                         0,
1435:                                         0,
1436:                                         0,
1437:                                 /* 94*/ 0,
1438:                                         0,
1439:                                         MatMatTransposeMultSymbolic_MPIDense_MPIDense,
1440:                                         MatMatTransposeMultNumeric_MPIDense_MPIDense,
1441:                                         0,
1442:                                 /* 99*/ MatProductSetFromOptions_MPIDense,
1443:                                         0,
1444:                                         0,
1445:                                         MatConjugate_MPIDense,
1446:                                         0,
1447:                                 /*104*/ 0,
1448:                                         MatRealPart_MPIDense,
1449:                                         MatImaginaryPart_MPIDense,
1450:                                         0,
1451:                                         0,
1452:                                 /*109*/ 0,
1453:                                         0,
1454:                                         0,
1455:                                         MatGetColumnVector_MPIDense,
1456:                                         MatMissingDiagonal_MPIDense,
1457:                                 /*114*/ 0,
1458:                                         0,
1459:                                         0,
1460:                                         0,
1461:                                         0,
1462:                                 /*119*/ 0,
1463:                                         0,
1464:                                         0,
1465:                                         0,
1466:                                         0,
1467:                                 /*124*/ 0,
1468:                                         MatGetColumnNorms_MPIDense,
1469:                                         0,
1470:                                         0,
1471:                                         0,
1472:                                 /*129*/ 0,
1473:                                         0,
1474:                                         MatTransposeMatMultSymbolic_MPIDense_MPIDense,
1475:                                         MatTransposeMatMultNumeric_MPIDense_MPIDense,
1476:                                         0,
1477:                                 /*134*/ 0,
1478:                                         0,
1479:                                         0,
1480:                                         0,
1481:                                         0,
1482:                                 /*139*/ 0,
1483:                                         0,
1484:                                         0,
1485:                                         0,
1486:                                         0,
1487:                                         MatCreateMPIMatConcatenateSeqMat_MPIDense,
1488:                                 /*145*/ 0,
1489:                                         0,
1490:                                         0
1491: };

1493: PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1494: {
1495:   Mat_MPIDense   *a = (Mat_MPIDense*)mat->data;
1496:   PetscBool      iscuda;

1500:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1501:   PetscLayoutSetUp(mat->rmap);
1502:   PetscLayoutSetUp(mat->cmap);
1503:   if (!a->A) {
1504:     MatCreate(PETSC_COMM_SELF,&a->A);
1505:     PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
1506:     MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);
1507:   }
1508:   PetscObjectTypeCompare((PetscObject)mat,MATMPIDENSECUDA,&iscuda);
1509:   MatSetType(a->A,iscuda ? MATSEQDENSECUDA : MATSEQDENSE);
1510:   MatSeqDenseSetPreallocation(a->A,data);
1511:   mat->preallocated = PETSC_TRUE;
1512:   return(0);
1513: }

1515: #if defined(PETSC_HAVE_ELEMENTAL)
1516: PETSC_INTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
1517: {
1518:   Mat            mat_elemental;
1520:   PetscScalar    *v;
1521:   PetscInt       m=A->rmap->n,N=A->cmap->N,rstart=A->rmap->rstart,i,*rows,*cols;

1524:   if (reuse == MAT_REUSE_MATRIX) {
1525:     mat_elemental = *newmat;
1526:     MatZeroEntries(*newmat);
1527:   } else {
1528:     MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);
1529:     MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);
1530:     MatSetType(mat_elemental,MATELEMENTAL);
1531:     MatSetUp(mat_elemental);
1532:     MatSetOption(mat_elemental,MAT_ROW_ORIENTED,PETSC_FALSE);
1533:   }

1535:   PetscMalloc2(m,&rows,N,&cols);
1536:   for (i=0; i<N; i++) cols[i] = i;
1537:   for (i=0; i<m; i++) rows[i] = rstart + i;

1539:   /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
1540:   MatDenseGetArray(A,&v);
1541:   MatSetValues(mat_elemental,m,rows,N,cols,v,ADD_VALUES);
1542:   MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);
1543:   MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);
1544:   MatDenseRestoreArray(A,&v);
1545:   PetscFree2(rows,cols);

1547:   if (reuse == MAT_INPLACE_MATRIX) {
1548:     MatHeaderReplace(A,&mat_elemental);
1549:   } else {
1550:     *newmat = mat_elemental;
1551:   }
1552:   return(0);
1553: }
1554: #endif

1556: static PetscErrorCode MatDenseGetColumn_MPIDense(Mat A,PetscInt col,PetscScalar **vals)
1557: {
1558:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;

1562:   MatDenseGetColumn(mat->A,col,vals);
1563:   return(0);
1564: }

1566: static PetscErrorCode MatDenseRestoreColumn_MPIDense(Mat A,PetscScalar **vals)
1567: {
1568:   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;

1572:   MatDenseRestoreColumn(mat->A,vals);
1573:   return(0);
1574: }

1576: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
1577: {
1579:   Mat_MPIDense   *mat;
1580:   PetscInt       m,nloc,N;

1583:   MatGetSize(inmat,&m,&N);
1584:   MatGetLocalSize(inmat,NULL,&nloc);
1585:   if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
1586:     PetscInt sum;

1588:     if (n == PETSC_DECIDE) {
1589:       PetscSplitOwnership(comm,&n,&N);
1590:     }
1591:     /* Check sum(n) = N */
1592:     MPIU_Allreduce(&n,&sum,1,MPIU_INT,MPI_SUM,comm);
1593:     if (sum != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Sum of local columns %D != global columns %D",sum,N);

1595:     MatCreateDense(comm,m,n,PETSC_DETERMINE,N,NULL,outmat);
1596:   }

1598:   /* numeric phase */
1599:   mat = (Mat_MPIDense*)(*outmat)->data;
1600:   MatCopy(inmat,mat->A,SAME_NONZERO_PATTERN);
1601:   MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY);
1602:   MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY);
1603:   return(0);
1604: }

1606: #if defined(PETSC_HAVE_CUDA)
1607: PetscErrorCode MatConvert_MPIDenseCUDA_MPIDense(Mat M,MatType type,MatReuse reuse,Mat *newmat)
1608: {
1609:   Mat            B;
1610:   Mat_MPIDense   *m;

1614:   if (reuse == MAT_INITIAL_MATRIX) {
1615:     MatDuplicate(M,MAT_COPY_VALUES,newmat);
1616:   } else if (reuse == MAT_REUSE_MATRIX) {
1617:     MatCopy(M,*newmat,SAME_NONZERO_PATTERN);
1618:   }

1620:   B    = *newmat;
1621:   MatBindToCPU_MPIDenseCUDA(B,PETSC_TRUE);
1622:   PetscFree(B->defaultvectype);
1623:   PetscStrallocpy(VECSTANDARD,&B->defaultvectype);
1624:   PetscObjectChangeTypeName((PetscObject)B,MATMPIDENSE);
1625:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpidensecuda_mpidense_C",NULL);
1626:   PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_mpiaij_mpidensecuda_C",NULL);
1627:   PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_mpiaijcusparse_mpidensecuda_C",NULL);
1628:   PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_mpidensecuda_mpiaij_C",NULL);
1629:   PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_mpidensecuda_mpiaijcusparse_C",NULL);
1630:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAGetArray_C",NULL);
1631:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAGetArrayRead_C",NULL);
1632:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAGetArrayWrite_C",NULL);
1633:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDARestoreArray_C",NULL);
1634:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDARestoreArrayRead_C",NULL);
1635:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDARestoreArrayWrite_C",NULL);
1636:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAPlaceArray_C",NULL);
1637:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAResetArray_C",NULL);
1638:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAReplaceArray_C",NULL);
1639:   m    = (Mat_MPIDense*)(B)->data;
1640:   if (m->A) {
1641:     MatConvert(m->A,MATSEQDENSE,MAT_INPLACE_MATRIX,&m->A);
1642:     MatSetUpMultiply_MPIDense(B);
1643:   }
1644:   B->ops->bindtocpu = NULL;
1645:   B->offloadmask    = PETSC_OFFLOAD_CPU;
1646:   return(0);
1647: }

1649: PetscErrorCode MatConvert_MPIDense_MPIDenseCUDA(Mat M,MatType type,MatReuse reuse,Mat *newmat)
1650: {
1651:   Mat            B;
1652:   Mat_MPIDense   *m;

1656:   if (reuse == MAT_INITIAL_MATRIX) {
1657:     MatDuplicate(M,MAT_COPY_VALUES,newmat);
1658:   } else if (reuse == MAT_REUSE_MATRIX) {
1659:     MatCopy(M,*newmat,SAME_NONZERO_PATTERN);
1660:   }

1662:   B    = *newmat;
1663:   PetscFree(B->defaultvectype);
1664:   PetscStrallocpy(VECCUDA,&B->defaultvectype);
1665:   PetscObjectChangeTypeName((PetscObject)B,MATMPIDENSECUDA);
1666:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpidensecuda_mpidense_C",                    MatConvert_MPIDenseCUDA_MPIDense);
1667:   PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_mpiaij_mpidensecuda_C",        MatProductSetFromOptions_MPIAIJ_MPIDense);
1668:   PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_mpiaijcusparse_mpidensecuda_C",MatProductSetFromOptions_MPIAIJ_MPIDense);
1669:   PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_mpidensecuda_mpiaij_C",        MatProductSetFromOptions_MPIDense_MPIAIJ);
1670:   PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_mpidensecuda_mpiaijcusparse_C",MatProductSetFromOptions_MPIDense_MPIAIJ);
1671:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAGetArray_C",                                MatDenseCUDAGetArray_MPIDenseCUDA);
1672:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAGetArrayRead_C",                            MatDenseCUDAGetArrayRead_MPIDenseCUDA);
1673:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAGetArrayWrite_C",                           MatDenseCUDAGetArrayWrite_MPIDenseCUDA);
1674:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDARestoreArray_C",                            MatDenseCUDARestoreArray_MPIDenseCUDA);
1675:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDARestoreArrayRead_C",                        MatDenseCUDARestoreArrayRead_MPIDenseCUDA);
1676:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDARestoreArrayWrite_C",                       MatDenseCUDARestoreArrayWrite_MPIDenseCUDA);
1677:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAPlaceArray_C",                              MatDenseCUDAPlaceArray_MPIDenseCUDA);
1678:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAResetArray_C",                              MatDenseCUDAResetArray_MPIDenseCUDA);
1679:   PetscObjectComposeFunction((PetscObject)B,"MatDenseCUDAReplaceArray_C",                            MatDenseCUDAReplaceArray_MPIDenseCUDA);
1680:   m    = (Mat_MPIDense*)(B)->data;
1681:   if (m->A) {
1682:     MatConvert(m->A,MATSEQDENSECUDA,MAT_INPLACE_MATRIX,&m->A);
1683:     MatSetUpMultiply_MPIDense(B);
1684:     B->offloadmask = PETSC_OFFLOAD_BOTH;
1685:   } else {
1686:     B->offloadmask = PETSC_OFFLOAD_UNALLOCATED;
1687:   }
1688:   MatBindToCPU_MPIDenseCUDA(B,PETSC_FALSE);

1690:   B->ops->bindtocpu = MatBindToCPU_MPIDenseCUDA;
1691:   return(0);
1692: }
1693: #endif

1695: PetscErrorCode MatDenseGetColumnVec_MPIDense(Mat A,PetscInt col,Vec *v)
1696: {
1697:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1699:   PetscInt       lda;

1702:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1703:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1704:   if (!a->cvec) {
1705:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,A->rmap->N,NULL,&a->cvec);
1706:   }
1707:   a->vecinuse = col + 1;
1708:   MatDenseGetLDA(a->A,&lda);
1709:   MatDenseGetArray(a->A,(PetscScalar**)&a->ptrinuse);
1710:   VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)lda);
1711:   *v   = a->cvec;
1712:   return(0);
1713: }

1715: PetscErrorCode MatDenseRestoreColumnVec_MPIDense(Mat A,PetscInt col,Vec *v)
1716: {
1717:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1721:   if (!a->vecinuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
1722:   if (!a->cvec) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column vector");
1723:   a->vecinuse = 0;
1724:   MatDenseRestoreArray(a->A,(PetscScalar**)&a->ptrinuse);
1725:   VecResetArray(a->cvec);
1726:   *v   = NULL;
1727:   return(0);
1728: }

1730: PetscErrorCode MatDenseGetColumnVecRead_MPIDense(Mat A,PetscInt col,Vec *v)
1731: {
1732:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1734:   PetscInt       lda;

1737:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1738:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1739:   if (!a->cvec) {
1740:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,A->rmap->N,NULL,&a->cvec);
1741:   }
1742:   a->vecinuse = col + 1;
1743:   MatDenseGetLDA(a->A,&lda);
1744:   MatDenseGetArrayRead(a->A,&a->ptrinuse);
1745:   VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)lda);
1746:   VecLockReadPush(a->cvec);
1747:   *v   = a->cvec;
1748:   return(0);
1749: }

1751: PetscErrorCode MatDenseRestoreColumnVecRead_MPIDense(Mat A,PetscInt col,Vec *v)
1752: {
1753:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1757:   if (!a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
1758:   if (!a->cvec) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing internal column vector");
1759:   a->vecinuse = 0;
1760:   MatDenseRestoreArrayRead(a->A,&a->ptrinuse);
1761:   VecLockReadPop(a->cvec);
1762:   VecResetArray(a->cvec);
1763:   *v   = NULL;
1764:   return(0);
1765: }

1767: PetscErrorCode MatDenseGetColumnVecWrite_MPIDense(Mat A,PetscInt col,Vec *v)
1768: {
1769:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1771:   PetscInt       lda;

1774:   if (a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1775:   if (a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1776:   if (!a->cvec) {
1777:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,A->rmap->N,NULL,&a->cvec);
1778:   }
1779:   a->vecinuse = col + 1;
1780:   MatDenseGetLDA(a->A,&lda);
1781:   MatDenseGetArrayWrite(a->A,(PetscScalar**)&a->ptrinuse);
1782:   VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)lda);
1783:   *v   = a->cvec;
1784:   return(0);
1785: }

1787: PetscErrorCode MatDenseRestoreColumnVecWrite_MPIDense(Mat A,PetscInt col,Vec *v)
1788: {
1789:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;

1793:   if (!a->vecinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first");
1794:   if (!a->cvec) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing internal column vector");
1795:   a->vecinuse = 0;
1796:   MatDenseRestoreArrayWrite(a->A,(PetscScalar**)&a->ptrinuse);
1797:   VecResetArray(a->cvec);
1798:   *v   = NULL;
1799:   return(0);
1800: }

1802: PetscErrorCode MatDenseGetSubMatrix_MPIDense(Mat A,PetscInt cbegin,PetscInt cend,Mat *v)
1803: {
1804:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1805:   Mat_MPIDense   *c;
1807:   MPI_Comm       comm;
1808:   PetscBool      setup = PETSC_FALSE;

1811:   PetscObjectGetComm((PetscObject)A,&comm);
1812:   if (a->vecinuse) SETERRQ(comm,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first");
1813:   if (a->matinuse) SETERRQ(comm,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1814:   if (!a->cmat) {
1815:     setup = PETSC_TRUE;
1816:     MatCreate(comm,&a->cmat);
1817:     PetscLogObjectParent((PetscObject)A,(PetscObject)a->cmat);
1818:     MatSetType(a->cmat,((PetscObject)A)->type_name);
1819:     PetscLayoutReference(A->rmap,&a->cmat->rmap);
1820:     PetscLayoutSetSize(a->cmat->cmap,cend-cbegin);
1821:     PetscLayoutSetUp(a->cmat->cmap);
1822:   } else if (cend-cbegin != a->cmat->cmap->N) {
1823:     setup = PETSC_TRUE;
1824:     PetscLayoutDestroy(&a->cmat->cmap);
1825:     PetscLayoutCreate(comm,&a->cmat->cmap);
1826:     PetscLayoutSetSize(a->cmat->cmap,cend-cbegin);
1827:     PetscLayoutSetUp(a->cmat->cmap);
1828:   }
1829:   c = (Mat_MPIDense*)a->cmat->data;
1830:   if (c->A) SETERRQ(comm,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first");
1831:   MatDenseGetSubMatrix(a->A,cbegin,cend,&c->A);
1832:   if (setup) { /* do we really need this? */
1833:     MatSetUpMultiply_MPIDense(a->cmat);
1834:   }
1835:   a->cmat->preallocated = PETSC_TRUE;
1836:   a->cmat->assembled = PETSC_TRUE;
1837:   a->matinuse = cbegin + 1;
1838:   *v = a->cmat;
1839:   return(0);
1840: }

1842: PetscErrorCode MatDenseRestoreSubMatrix_MPIDense(Mat A,Mat *v)
1843: {
1844:   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1845:   Mat_MPIDense   *c;

1849:   if (!a->matinuse) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Need to call MatDenseGetSubMatrix() first");
1850:   if (!a->cmat) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing internal matrix");
1851:   if (*v != a->cmat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not the matrix obtained from MatDenseGetSubMatrix()");
1852:   a->matinuse = 0;
1853:   c    = (Mat_MPIDense*)a->cmat->data;
1854:   MatDenseRestoreSubMatrix(a->A,&c->A);
1855:   *v   = NULL;
1856:   return(0);
1857: }

1859: PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat)
1860: {
1861:   Mat_MPIDense   *a;

1865:   PetscNewLog(mat,&a);
1866:   mat->data = (void*)a;
1867:   PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));

1869:   mat->insertmode = NOT_SET_VALUES;

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

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

1876:   /* stuff used for matrix vector multiply */
1877:   a->lvec        = 0;
1878:   a->Mvctx       = 0;
1879:   a->roworiented = PETSC_TRUE;

1881:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",MatDenseGetLDA_MPIDense);
1882:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseSetLDA_C",MatDenseSetLDA_MPIDense);
1883:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);
1884:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);
1885:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",MatDenseGetArrayRead_MPIDense);
1886:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",MatDenseRestoreArrayRead_MPIDense);
1887:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayWrite_C",MatDenseGetArrayWrite_MPIDense);
1888:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayWrite_C",MatDenseRestoreArrayWrite_MPIDense);
1889:   PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",MatDensePlaceArray_MPIDense);
1890:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",MatDenseResetArray_MPIDense);
1891:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseReplaceArray_C",MatDenseReplaceArray_MPIDense);
1892:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVec_C",MatDenseGetColumnVec_MPIDense);
1893:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVec_C",MatDenseRestoreColumnVec_MPIDense);
1894:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecRead_C",MatDenseGetColumnVecRead_MPIDense);
1895:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecRead_C",MatDenseRestoreColumnVecRead_MPIDense);
1896:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecWrite_C",MatDenseGetColumnVecWrite_MPIDense);
1897:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecWrite_C",MatDenseRestoreColumnVecWrite_MPIDense);
1898:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetSubMatrix_C",MatDenseGetSubMatrix_MPIDense);
1899:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreSubMatrix_C",MatDenseRestoreSubMatrix_MPIDense);
1900: #if defined(PETSC_HAVE_ELEMENTAL)
1901:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",MatConvert_MPIDense_Elemental);
1902: #endif
1903: #if defined(PETSC_HAVE_SCALAPACK)
1904:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_scalapack_C",MatConvert_Dense_ScaLAPACK);
1905: #endif
1906: #if defined(PETSC_HAVE_CUDA)
1907:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_mpidensecuda_C",MatConvert_MPIDense_MPIDenseCUDA);
1908: #endif
1909:   PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);
1910:   PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_mpiaij_mpidense_C",MatProductSetFromOptions_MPIAIJ_MPIDense);
1911:   PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_mpidense_mpiaij_C",MatProductSetFromOptions_MPIDense_MPIAIJ);
1912: #if defined(PETSC_HAVE_CUDA)
1913:   PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_mpiaijcusparse_mpidense_C",MatProductSetFromOptions_MPIAIJ_MPIDense);
1914:   PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_mpidense_mpiaijcusparse_C",MatProductSetFromOptions_MPIDense_MPIAIJ);
1915: #endif

1917:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",MatDenseGetColumn_MPIDense);
1918:   PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",MatDenseRestoreColumn_MPIDense);
1919:   PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);
1920:   return(0);
1921: }

1923: /*MC
1924:    MATMPIDENSECUDA - MATMPIDENSECUDA = "mpidensecuda" - A matrix type to be used for distributed dense matrices on GPUs.

1926:    Options Database Keys:
1927: . -mat_type mpidensecuda - sets the matrix type to "mpidensecuda" during a call to MatSetFromOptions()

1929:   Level: beginner

1931: .seealso:

1933: M*/
1934: #if defined(PETSC_HAVE_CUDA)
1935: PETSC_EXTERN PetscErrorCode MatCreate_MPIDenseCUDA(Mat B)
1936: {

1940:   MatCreate_MPIDense(B);
1941:   MatConvert_MPIDense_MPIDenseCUDA(B,MATMPIDENSECUDA,MAT_INPLACE_MATRIX,&B);
1942:   return(0);
1943: }
1944: #endif

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

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

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

1955:   Level: beginner


1958: .seealso: MATSEQDENSE,MATMPIDENSE,MATDENSECUDA
1959: M*/

1961: /*MC
1962:    MATDENSECUDA - MATDENSECUDA = "densecuda" - A matrix type to be used for dense matrices on GPUs.

1964:    This matrix type is identical to MATSEQDENSECUDA when constructed with a single process communicator,
1965:    and MATMPIDENSECUDA otherwise.

1967:    Options Database Keys:
1968: . -mat_type densecuda - sets the matrix type to "densecuda" during a call to MatSetFromOptions()

1970:   Level: beginner

1972: .seealso: MATSEQDENSECUDA,MATMPIDENSECUDA,MATDENSE
1973: M*/

1975: /*@C
1976:    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries

1978:    Collective

1980:    Input Parameters:
1981: .  B - the matrix
1982: -  data - optional location of matrix data.  Set data=NULL for PETSc
1983:    to control all matrix memory allocation.

1985:    Notes:
1986:    The dense format is fully compatible with standard Fortran 77
1987:    storage by columns.

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

1993:    Level: intermediate

1995: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1996: @*/
1997: PetscErrorCode  MatMPIDenseSetPreallocation(Mat B,PetscScalar *data)
1998: {

2003:   PetscTryMethod(B,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(B,data));
2004:   return(0);
2005: }

2007: /*@
2008:    MatDensePlaceArray - Allows one to replace the array in a dense matrix with an
2009:    array provided by the user. This is useful to avoid copying an array
2010:    into a matrix

2012:    Not Collective

2014:    Input Parameters:
2015: +  mat - the matrix
2016: -  array - the array in column major order

2018:    Notes:
2019:    You can return to the original array with a call to MatDenseResetArray(). The user is responsible for freeing this array; it will not be
2020:    freed when the matrix is destroyed.

2022:    Level: developer

2024: .seealso: MatDenseGetArray(), MatDenseResetArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()

2026: @*/
2027: PetscErrorCode  MatDensePlaceArray(Mat mat,const PetscScalar *array)
2028: {

2033:   PetscUseMethod(mat,"MatDensePlaceArray_C",(Mat,const PetscScalar*),(mat,array));
2034:   PetscObjectStateIncrease((PetscObject)mat);
2035: #if defined(PETSC_HAVE_CUDA)
2036:   mat->offloadmask = PETSC_OFFLOAD_CPU;
2037: #endif
2038:   return(0);
2039: }

2041: /*@
2042:    MatDenseResetArray - Resets the matrix array to that it previously had before the call to MatDensePlaceArray()

2044:    Not Collective

2046:    Input Parameters:
2047: .  mat - the matrix

2049:    Notes:
2050:    You can only call this after a call to MatDensePlaceArray()

2052:    Level: developer

2054: .seealso: MatDenseGetArray(), MatDensePlaceArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()

2056: @*/
2057: PetscErrorCode  MatDenseResetArray(Mat mat)
2058: {

2063:   PetscUseMethod(mat,"MatDenseResetArray_C",(Mat),(mat));
2064:   PetscObjectStateIncrease((PetscObject)mat);
2065:   return(0);
2066: }

2068: /*@
2069:    MatDenseReplaceArray - Allows one to replace the array in a dense matrix with an
2070:    array provided by the user. This is useful to avoid copying an array
2071:    into a matrix

2073:    Not Collective

2075:    Input Parameters:
2076: +  mat - the matrix
2077: -  array - the array in column major order

2079:    Notes:
2080:    The memory passed in MUST be obtained with PetscMalloc() and CANNOT be
2081:    freed by the user. It will be freed when the matrix is destroyed.

2083:    Level: developer

2085: .seealso: MatDenseGetArray(), VecReplaceArray()
2086: @*/
2087: PetscErrorCode  MatDenseReplaceArray(Mat mat,const PetscScalar *array)
2088: {

2093:   PetscUseMethod(mat,"MatDenseReplaceArray_C",(Mat,const PetscScalar*),(mat,array));
2094:   PetscObjectStateIncrease((PetscObject)mat);
2095: #if defined(PETSC_HAVE_CUDA)
2096:   mat->offloadmask = PETSC_OFFLOAD_CPU;
2097: #endif
2098:   return(0);
2099: }

2101: #if defined(PETSC_HAVE_CUDA)
2102: /*@C
2103:    MatDenseCUDAPlaceArray - Allows one to replace the GPU array in a dense matrix with an
2104:    array provided by the user. This is useful to avoid copying an array
2105:    into a matrix

2107:    Not Collective

2109:    Input Parameters:
2110: +  mat - the matrix
2111: -  array - the array in column major order

2113:    Notes:
2114:    You can return to the original array with a call to MatDenseCUDAResetArray(). The user is responsible for freeing this array; it will not be
2115:    freed when the matrix is destroyed. The array must have been allocated with cudaMalloc().

2117:    Level: developer

2119: .seealso: MatDenseCUDAGetArray(), MatDenseCUDAResetArray()
2120: @*/
2121: PetscErrorCode  MatDenseCUDAPlaceArray(Mat mat,const PetscScalar *array)
2122: {

2127:   PetscUseMethod(mat,"MatDenseCUDAPlaceArray_C",(Mat,const PetscScalar*),(mat,array));
2128:   PetscObjectStateIncrease((PetscObject)mat);
2129:   mat->offloadmask = PETSC_OFFLOAD_GPU;
2130:   return(0);
2131: }

2133: /*@C
2134:    MatDenseCUDAResetArray - Resets the matrix array to that it previously had before the call to MatDenseCUDAPlaceArray()

2136:    Not Collective

2138:    Input Parameters:
2139: .  mat - the matrix

2141:    Notes:
2142:    You can only call this after a call to MatDenseCUDAPlaceArray()

2144:    Level: developer

2146: .seealso: MatDenseCUDAGetArray(), MatDenseCUDAPlaceArray()

2148: @*/
2149: PetscErrorCode  MatDenseCUDAResetArray(Mat mat)
2150: {

2155:   PetscUseMethod(mat,"MatDenseCUDAResetArray_C",(Mat),(mat));
2156:   PetscObjectStateIncrease((PetscObject)mat);
2157:   return(0);
2158: }

2160: /*@C
2161:    MatDenseCUDAReplaceArray - Allows one to replace the GPU array in a dense matrix with an
2162:    array provided by the user. This is useful to avoid copying an array
2163:    into a matrix

2165:    Not Collective

2167:    Input Parameters:
2168: +  mat - the matrix
2169: -  array - the array in column major order

2171:    Notes:
2172:    This permanently replaces the GPU array and frees the memory associated with the old GPU array.
2173:    The memory passed in CANNOT be freed by the user. It will be freed
2174:    when the matrix is destroyed. The array should respect the matrix leading dimension.

2176:    Level: developer

2178: .seealso: MatDenseCUDAGetArray(), MatDenseCUDAPlaceArray(), MatDenseCUDAResetArray()
2179: @*/
2180: PetscErrorCode  MatDenseCUDAReplaceArray(Mat mat,const PetscScalar *array)
2181: {

2186:   PetscUseMethod(mat,"MatDenseCUDAReplaceArray_C",(Mat,const PetscScalar*),(mat,array));
2187:   PetscObjectStateIncrease((PetscObject)mat);
2188:   mat->offloadmask = PETSC_OFFLOAD_GPU;
2189:   return(0);
2190: }

2192: /*@C
2193:    MatDenseCUDAGetArrayWrite - Provides write access to the CUDA buffer inside a dense matrix.

2195:    Not Collective

2197:    Input Parameters:
2198: .  A - the matrix

2200:    Output Parameters
2201: .  array - the GPU array in column major order

2203:    Notes:
2204:    The data on the GPU may not be updated due to operations done on the CPU. If you need updated data, use MatDenseCUDAGetArray(). The array must be restored with MatDenseCUDARestoreArrayWrite() when no longer needed.

2206:    Level: developer

2208: .seealso: MatDenseCUDAGetArray(), MatDenseCUDARestoreArray(), MatDenseCUDARestoreArrayWrite(), MatDenseCUDAGetArrayRead(), MatDenseCUDARestoreArrayRead()
2209: @*/
2210: PetscErrorCode MatDenseCUDAGetArrayWrite(Mat A, PetscScalar **a)
2211: {

2216:   PetscUseMethod(A,"MatDenseCUDAGetArrayWrite_C",(Mat,PetscScalar**),(A,a));
2217:   PetscObjectStateIncrease((PetscObject)A);
2218:   return(0);
2219: }

2221: /*@C
2222:    MatDenseCUDARestoreArrayWrite - Restore write access to the CUDA buffer inside a dense matrix previously obtained with MatDenseCUDAGetArrayWrite().

2224:    Not Collective

2226:    Input Parameters:
2227: +  A - the matrix
2228: -  array - the GPU array in column major order

2230:    Notes:

2232:    Level: developer

2234: .seealso: MatDenseCUDAGetArray(), MatDenseCUDARestoreArray(), MatDenseCUDAGetArrayWrite(), MatDenseCUDARestoreArrayRead(), MatDenseCUDAGetArrayRead()
2235: @*/
2236: PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat A, PetscScalar **a)
2237: {

2242:   PetscUseMethod(A,"MatDenseCUDARestoreArrayWrite_C",(Mat,PetscScalar**),(A,a));
2243:   PetscObjectStateIncrease((PetscObject)A);
2244:   A->offloadmask = PETSC_OFFLOAD_GPU;
2245:   return(0);
2246: }

2248: /*@C
2249:    MatDenseCUDAGetArrayRead - Provides read-only access to the CUDA buffer inside a dense matrix. The array must be restored with MatDenseCUDARestoreArrayRead() when no longer needed.

2251:    Not Collective

2253:    Input Parameters:
2254: .  A - the matrix

2256:    Output Parameters
2257: .  array - the GPU array in column major order

2259:    Notes:
2260:    Data can be copied to the GPU due to operations done on the CPU. If you need write only access, use MatDenseCUDAGetArrayWrite().

2262:    Level: developer

2264: .seealso: MatDenseCUDAGetArray(), MatDenseCUDARestoreArray(), MatDenseCUDARestoreArrayWrite(), MatDenseCUDAGetArrayWrite(), MatDenseCUDARestoreArrayRead()
2265: @*/
2266: PetscErrorCode MatDenseCUDAGetArrayRead(Mat A, const PetscScalar **a)
2267: {

2272:   PetscUseMethod(A,"MatDenseCUDAGetArrayRead_C",(Mat,const PetscScalar**),(A,a));
2273:   return(0);
2274: }

2276: /*@C
2277:    MatDenseCUDARestoreArrayRead - Restore read-only access to the CUDA buffer inside a dense matrix previously obtained with a call to MatDenseCUDAGetArrayRead().

2279:    Not Collective

2281:    Input Parameters:
2282: +  A - the matrix
2283: -  array - the GPU array in column major order

2285:    Notes:
2286:    Data can be copied to the GPU due to operations done on the CPU. If you need write only access, use MatDenseCUDAGetArrayWrite().

2288:    Level: developer

2290: .seealso: MatDenseCUDAGetArray(), MatDenseCUDARestoreArray(), MatDenseCUDARestoreArrayWrite(), MatDenseCUDAGetArrayWrite(), MatDenseCUDAGetArrayRead()
2291: @*/
2292: PetscErrorCode MatDenseCUDARestoreArrayRead(Mat A, const PetscScalar **a)
2293: {

2297:   PetscUseMethod(A,"MatDenseCUDARestoreArrayRead_C",(Mat,const PetscScalar**),(A,a));
2298:   return(0);
2299: }

2301: /*@C
2302:    MatDenseCUDAGetArray - Provides access to the CUDA buffer inside a dense matrix. The array must be restored with MatDenseCUDARestoreArray() when no longer needed.

2304:    Not Collective

2306:    Input Parameters:
2307: .  A - the matrix

2309:    Output Parameters
2310: .  array - the GPU array in column major order

2312:    Notes:
2313:    Data can be copied to the GPU due to operations done on the CPU. If you need write only access, use MatDenseCUDAGetArrayWrite(). For read-only access, use MatDenseCUDAGetArrayRead().

2315:    Level: developer

2317: .seealso: MatDenseCUDAGetArrayRead(), MatDenseCUDARestoreArray(), MatDenseCUDARestoreArrayWrite(), MatDenseCUDAGetArrayWrite(), MatDenseCUDARestoreArrayRead()
2318: @*/
2319: PetscErrorCode MatDenseCUDAGetArray(Mat A, PetscScalar **a)
2320: {

2325:   PetscUseMethod(A,"MatDenseCUDAGetArray_C",(Mat,PetscScalar**),(A,a));
2326:   PetscObjectStateIncrease((PetscObject)A);
2327:   return(0);
2328: }

2330: /*@C
2331:    MatDenseCUDARestoreArray - Restore access to the CUDA buffer inside a dense matrix previously obtained with MatDenseCUDAGetArray().

2333:    Not Collective

2335:    Input Parameters:
2336: +  A - the matrix
2337: -  array - the GPU array in column major order

2339:    Notes:

2341:    Level: developer

2343: .seealso: MatDenseCUDAGetArray(), MatDenseCUDARestoreArrayWrite(), MatDenseCUDAGetArrayWrite(), MatDenseCUDARestoreArrayRead(), MatDenseCUDAGetArrayRead()
2344: @*/
2345: PetscErrorCode MatDenseCUDARestoreArray(Mat A, PetscScalar **a)
2346: {

2351:   PetscUseMethod(A,"MatDenseCUDARestoreArray_C",(Mat,PetscScalar**),(A,a));
2352:   PetscObjectStateIncrease((PetscObject)A);
2353:   A->offloadmask = PETSC_OFFLOAD_GPU;
2354:   return(0);
2355: }
2356: #endif

2358: /*@C
2359:    MatCreateDense - Creates a matrix in dense format.

2361:    Collective

2363:    Input Parameters:
2364: +  comm - MPI communicator
2365: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2366: .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2367: .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
2368: .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
2369: -  data - optional location of matrix data.  Set data=NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
2370:    to control all matrix memory allocation.

2372:    Output Parameter:
2373: .  A - the matrix

2375:    Notes:
2376:    The dense format is fully compatible with standard Fortran 77
2377:    storage by columns.

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

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

2386:    Level: intermediate

2388: .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
2389: @*/
2390: PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
2391: {
2393:   PetscMPIInt    size;

2396:   MatCreate(comm,A);
2397:   MatSetSizes(*A,m,n,M,N);
2398:   MPI_Comm_size(comm,&size);
2399:   if (size > 1) {
2400:     PetscBool havedata = (PetscBool)!!data;

2402:     MatSetType(*A,MATMPIDENSE);
2403:     MatMPIDenseSetPreallocation(*A,data);
2404:     MPIU_Allreduce(MPI_IN_PLACE,&havedata,1,MPIU_BOOL,MPI_LOR,comm);
2405:     if (havedata) {  /* user provided data array, so no need to assemble */
2406:       MatSetUpMultiply_MPIDense(*A);
2407:       (*A)->assembled = PETSC_TRUE;
2408:     }
2409:   } else {
2410:     MatSetType(*A,MATSEQDENSE);
2411:     MatSeqDenseSetPreallocation(*A,data);
2412:   }
2413:   return(0);
2414: }

2416: #if defined(PETSC_HAVE_CUDA)
2417: /*@C
2418:    MatCreateDenseCUDA - Creates a matrix in dense format using CUDA.

2420:    Collective

2422:    Input Parameters:
2423: +  comm - MPI communicator
2424: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2425: .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2426: .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
2427: .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
2428: -  data - optional location of GPU matrix data.  Set data=NULL for PETSc
2429:    to control matrix memory allocation.

2431:    Output Parameter:
2432: .  A - the matrix

2434:    Notes:

2436:    Level: intermediate

2438: .seealso: MatCreate(), MatCreateDense()
2439: @*/
2440: PetscErrorCode  MatCreateDenseCUDA(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
2441: {
2443:   PetscMPIInt    size;

2446:   MatCreate(comm,A);
2448:   MatSetSizes(*A,m,n,M,N);
2449:   MPI_Comm_size(comm,&size);
2450:   if (size > 1) {
2451:     MatSetType(*A,MATMPIDENSECUDA);
2452:     MatMPIDenseCUDASetPreallocation(*A,data);
2453:     if (data) {  /* user provided data array, so no need to assemble */
2454:       MatSetUpMultiply_MPIDense(*A);
2455:       (*A)->assembled = PETSC_TRUE;
2456:     }
2457:   } else {
2458:     MatSetType(*A,MATSEQDENSECUDA);
2459:     MatSeqDenseCUDASetPreallocation(*A,data);
2460:   }
2461:   return(0);
2462: }
2463: #endif

2465: static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
2466: {
2467:   Mat            mat;
2468:   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;

2472:   *newmat = 0;
2473:   MatCreate(PetscObjectComm((PetscObject)A),&mat);
2474:   MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2475:   MatSetType(mat,((PetscObject)A)->type_name);
2476:   a       = (Mat_MPIDense*)mat->data;

2478:   mat->factortype   = A->factortype;
2479:   mat->assembled    = PETSC_TRUE;
2480:   mat->preallocated = PETSC_TRUE;

2482:   mat->insertmode = NOT_SET_VALUES;
2483:   a->donotstash   = oldmat->donotstash;

2485:   PetscLayoutReference(A->rmap,&mat->rmap);
2486:   PetscLayoutReference(A->cmap,&mat->cmap);

2488:   MatDuplicate(oldmat->A,cpvalues,&a->A);
2489:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
2490:   MatSetUpMultiply_MPIDense(mat);

2492:   *newmat = mat;
2493:   return(0);
2494: }

2496: PetscErrorCode MatLoad_MPIDense(Mat newMat, PetscViewer viewer)
2497: {
2499:   PetscBool      isbinary;
2500: #if defined(PETSC_HAVE_HDF5)
2501:   PetscBool      ishdf5;
2502: #endif

2507:   /* force binary viewer to load .info file if it has not yet done so */
2508:   PetscViewerSetUp(viewer);
2509:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
2510: #if defined(PETSC_HAVE_HDF5)
2511:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,  &ishdf5);
2512: #endif
2513:   if (isbinary) {
2514:     MatLoad_Dense_Binary(newMat,viewer);
2515: #if defined(PETSC_HAVE_HDF5)
2516:   } else if (ishdf5) {
2517:     MatLoad_Dense_HDF5(newMat,viewer);
2518: #endif
2519:   } else SETERRQ2(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newMat)->type_name);
2520:   return(0);
2521: }

2523: static PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
2524: {
2525:   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
2526:   Mat            a,b;
2527:   PetscBool      flg;

2531:   a    = matA->A;
2532:   b    = matB->A;
2533:   MatEqual(a,b,&flg);
2534:   MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
2535:   return(0);
2536: }

2538: PetscErrorCode MatDestroy_MatTransMatMult_MPIDense_MPIDense(void *data)
2539: {
2540:   PetscErrorCode        ierr;
2541:   Mat_TransMatMultDense *atb = (Mat_TransMatMultDense *)data;

2544:   PetscFree2(atb->sendbuf,atb->recvcounts);
2545:   MatDestroy(&atb->atb);
2546:   PetscFree(atb);
2547:   return(0);
2548: }

2550: PetscErrorCode MatDestroy_MatMatTransMult_MPIDense_MPIDense(void *data)
2551: {
2552:   PetscErrorCode        ierr;
2553:   Mat_MatTransMultDense *abt = (Mat_MatTransMultDense *)data;

2556:   PetscFree2(abt->buf[0],abt->buf[1]);
2557:   PetscFree2(abt->recvcounts,abt->recvdispls);
2558:   PetscFree(abt);
2559:   return(0);
2560: }

2562: static PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C)
2563: {
2564:   Mat_MPIDense          *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data;
2565:   Mat_TransMatMultDense *atb;
2566:   PetscErrorCode        ierr;
2567:   MPI_Comm              comm;
2568:   PetscMPIInt           size,*recvcounts;
2569:   PetscScalar           *carray,*sendbuf;
2570:   const PetscScalar     *atbarray;
2571:   PetscInt              i,cN=C->cmap->N,cM=C->rmap->N,proc,k,j;
2572:   const PetscInt        *ranges;

2575:   MatCheckProduct(C,3);
2576:   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
2577:   atb = (Mat_TransMatMultDense *)C->product->data;
2578:   recvcounts = atb->recvcounts;
2579:   sendbuf = atb->sendbuf;

2581:   PetscObjectGetComm((PetscObject)A,&comm);
2582:   MPI_Comm_size(comm,&size);

2584:   /* compute atbarray = aseq^T * bseq */
2585:   MatTransposeMatMult(a->A,b->A,atb->atb ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX,PETSC_DEFAULT,&atb->atb);

2587:   MatGetOwnershipRanges(C,&ranges);

2589:   /* arrange atbarray into sendbuf */
2590:   MatDenseGetArrayRead(atb->atb,&atbarray);
2591:   for (proc=0, k=0; proc<size; proc++) {
2592:     for (j=0; j<cN; j++) {
2593:       for (i=ranges[proc]; i<ranges[proc+1]; i++) sendbuf[k++] = atbarray[i+j*cM];
2594:     }
2595:   }
2596:   MatDenseRestoreArrayRead(atb->atb,&atbarray);

2598:   /* sum all atbarray to local values of C */
2599:   MatDenseGetArrayWrite(c->A,&carray);
2600:   MPI_Reduce_scatter(sendbuf,carray,recvcounts,MPIU_SCALAR,MPIU_SUM,comm);
2601:   MatDenseRestoreArrayWrite(c->A,&carray);
2602:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2603:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2604:   return(0);
2605: }

2607: static PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
2608: {
2609:   PetscErrorCode        ierr;
2610:   MPI_Comm              comm;
2611:   PetscMPIInt           size;
2612:   PetscInt              cm=A->cmap->n,cM,cN=B->cmap->N;
2613:   Mat_TransMatMultDense *atb;
2614:   PetscBool             cisdense;
2615:   PetscInt              i;
2616:   const PetscInt        *ranges;

2619:   MatCheckProduct(C,3);
2620:   if (C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
2621:   PetscObjectGetComm((PetscObject)A,&comm);
2622:   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) {
2623:     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);
2624:   }

2626:   /* create matrix product C */
2627:   MatSetSizes(C,cm,B->cmap->n,A->cmap->N,B->cmap->N);
2628:   PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATMPIDENSE,MATMPIDENSECUDA,"");
2629:   if (!cisdense) {
2630:     MatSetType(C,((PetscObject)A)->type_name);
2631:   }
2632:   MatSetUp(C);

2634:   /* create data structure for reuse C */
2635:   MPI_Comm_size(comm,&size);
2636:   PetscNew(&atb);
2637:   cM   = C->rmap->N;
2638:   PetscMalloc2((size_t)cM*(size_t)cN,&atb->sendbuf,size,&atb->recvcounts);
2639:   MatGetOwnershipRanges(C,&ranges);
2640:   for (i=0; i<size; i++) atb->recvcounts[i] = (ranges[i+1] - ranges[i])*cN;

2642:   C->product->data    = atb;
2643:   C->product->destroy = MatDestroy_MatTransMatMult_MPIDense_MPIDense;
2644:   return(0);
2645: }

2647: static PetscErrorCode MatMatTransposeMultSymbolic_MPIDense_MPIDense(Mat A, Mat B, PetscReal fill, Mat C)
2648: {
2649:   PetscErrorCode        ierr;
2650:   MPI_Comm              comm;
2651:   PetscMPIInt           i, size;
2652:   PetscInt              maxRows, bufsiz;
2653:   PetscMPIInt           tag;
2654:   PetscInt              alg;
2655:   Mat_MatTransMultDense *abt;
2656:   Mat_Product           *product = C->product;
2657:   PetscBool             flg;

2660:   MatCheckProduct(C,4);
2661:   if (C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
2662:   /* check local size of A and B */
2663:   if (A->cmap->n != B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local column dimensions are incompatible, A (%D) != B (%D)",A->cmap->n,B->cmap->n);

2665:   PetscStrcmp(product->alg,"allgatherv",&flg);
2666:   alg  = flg ? 0 : 1;

2668:   /* setup matrix product C */
2669:   MatSetSizes(C,A->rmap->n,B->rmap->n,A->rmap->N,B->rmap->N);
2670:   MatSetType(C,MATMPIDENSE);
2671:   MatSetUp(C);
2672:   PetscObjectGetNewTag((PetscObject)C,&tag);

2674:   /* create data structure for reuse C */
2675:   PetscObjectGetComm((PetscObject)C,&comm);
2676:   MPI_Comm_size(comm,&size);
2677:   PetscNew(&abt);
2678:   abt->tag = tag;
2679:   abt->alg = alg;
2680:   switch (alg) {
2681:   case 1: /* alg: "cyclic" */
2682:     for (maxRows = 0, i = 0; i < size; i++) maxRows = PetscMax(maxRows, (B->rmap->range[i + 1] - B->rmap->range[i]));
2683:     bufsiz = A->cmap->N * maxRows;
2684:     PetscMalloc2(bufsiz,&(abt->buf[0]),bufsiz,&(abt->buf[1]));
2685:     break;
2686:   default: /* alg: "allgatherv" */
2687:     PetscMalloc2(B->rmap->n * B->cmap->N, &(abt->buf[0]), B->rmap->N * B->cmap->N, &(abt->buf[1]));
2688:     PetscMalloc2(size,&(abt->recvcounts),size+1,&(abt->recvdispls));
2689:     for (i = 0; i <= size; i++) abt->recvdispls[i] = B->rmap->range[i] * A->cmap->N;
2690:     for (i = 0; i < size; i++) abt->recvcounts[i] = abt->recvdispls[i + 1] - abt->recvdispls[i];
2691:     break;
2692:   }

2694:   C->product->data    = abt;
2695:   C->product->destroy = MatDestroy_MatMatTransMult_MPIDense_MPIDense;
2696:   C->ops->mattransposemultnumeric = MatMatTransposeMultNumeric_MPIDense_MPIDense;
2697:   return(0);
2698: }

2700: static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense_Cyclic(Mat A, Mat B, Mat C)
2701: {
2702:   Mat_MPIDense          *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data;
2703:   Mat_MatTransMultDense *abt;
2704:   PetscErrorCode        ierr;
2705:   MPI_Comm              comm;
2706:   PetscMPIInt           rank,size, sendsiz, recvsiz, sendto, recvfrom, recvisfrom;
2707:   PetscScalar           *sendbuf, *recvbuf=0, *cv;
2708:   PetscInt              i,cK=A->cmap->N,k,j,bn;
2709:   PetscScalar           _DOne=1.0,_DZero=0.0;
2710:   const PetscScalar     *av,*bv;
2711:   PetscBLASInt          cm, cn, ck, alda, blda, clda;
2712:   MPI_Request           reqs[2];
2713:   const PetscInt        *ranges;

2716:   MatCheckProduct(C,3);
2717:   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
2718:   abt  = (Mat_MatTransMultDense*)C->product->data;
2719:   PetscObjectGetComm((PetscObject)C,&comm);
2720:   MPI_Comm_rank(comm,&rank);
2721:   MPI_Comm_size(comm,&size);
2722:   MatDenseGetArrayRead(a->A,&av);
2723:   MatDenseGetArrayRead(b->A,&bv);
2724:   MatDenseGetArrayWrite(c->A,&cv);
2725:   MatDenseGetLDA(a->A,&i);
2726:   PetscBLASIntCast(i,&alda);
2727:   MatDenseGetLDA(b->A,&i);
2728:   PetscBLASIntCast(i,&blda);
2729:   MatDenseGetLDA(c->A,&i);
2730:   PetscBLASIntCast(i,&clda);
2731:   MatGetOwnershipRanges(B,&ranges);
2732:   bn   = B->rmap->n;
2733:   if (blda == bn) {
2734:     sendbuf = (PetscScalar*)bv;
2735:   } else {
2736:     sendbuf = abt->buf[0];
2737:     for (k = 0, i = 0; i < cK; i++) {
2738:       for (j = 0; j < bn; j++, k++) {
2739:         sendbuf[k] = bv[i * blda + j];
2740:       }
2741:     }
2742:   }
2743:   if (size > 1) {
2744:     sendto = (rank + size - 1) % size;
2745:     recvfrom = (rank + size + 1) % size;
2746:   } else {
2747:     sendto = recvfrom = 0;
2748:   }
2749:   PetscBLASIntCast(cK,&ck);
2750:   PetscBLASIntCast(c->A->rmap->n,&cm);
2751:   recvisfrom = rank;
2752:   for (i = 0; i < size; i++) {
2753:     /* we have finished receiving in sending, bufs can be read/modified */
2754:     PetscInt nextrecvisfrom = (recvisfrom + 1) % size; /* which process the next recvbuf will originate on */
2755:     PetscInt nextbn = ranges[nextrecvisfrom + 1] - ranges[nextrecvisfrom];

2757:     if (nextrecvisfrom != rank) {
2758:       /* start the cyclic sends from sendbuf, to recvbuf (which will switch to sendbuf) */
2759:       sendsiz = cK * bn;
2760:       recvsiz = cK * nextbn;
2761:       recvbuf = (i & 1) ? abt->buf[0] : abt->buf[1];
2762:       MPI_Isend(sendbuf, sendsiz, MPIU_SCALAR, sendto, abt->tag, comm, &reqs[0]);
2763:       MPI_Irecv(recvbuf, recvsiz, MPIU_SCALAR, recvfrom, abt->tag, comm, &reqs[1]);
2764:     }

2766:     /* local aseq * sendbuf^T */
2767:     PetscBLASIntCast(ranges[recvisfrom + 1] - ranges[recvisfrom], &cn);
2768:     if (cm && cn && ck) PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&cm,&cn,&ck,&_DOne,av,&alda,sendbuf,&cn,&_DZero,cv + clda * ranges[recvisfrom],&clda));

2770:     if (nextrecvisfrom != rank) {
2771:       /* wait for the sends and receives to complete, swap sendbuf and recvbuf */
2772:       MPI_Waitall(2, reqs, MPI_STATUSES_IGNORE);
2773:     }
2774:     bn = nextbn;
2775:     recvisfrom = nextrecvisfrom;
2776:     sendbuf = recvbuf;
2777:   }
2778:   MatDenseRestoreArrayRead(a->A,&av);
2779:   MatDenseRestoreArrayRead(b->A,&bv);
2780:   MatDenseRestoreArrayWrite(c->A,&cv);
2781:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2782:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2783:   return(0);
2784: }

2786: static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense_Allgatherv(Mat A, Mat B, Mat C)
2787: {
2788:   Mat_MPIDense          *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data;
2789:   Mat_MatTransMultDense *abt;
2790:   PetscErrorCode        ierr;
2791:   MPI_Comm              comm;
2792:   PetscMPIInt           size;
2793:   PetscScalar           *cv, *sendbuf, *recvbuf;
2794:   const PetscScalar     *av,*bv;
2795:   PetscInt              blda,i,cK=A->cmap->N,k,j,bn;
2796:   PetscScalar           _DOne=1.0,_DZero=0.0;
2797:   PetscBLASInt          cm, cn, ck, alda, clda;

2800:   MatCheckProduct(C,3);
2801:   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
2802:   abt  = (Mat_MatTransMultDense*)C->product->data;
2803:   PetscObjectGetComm((PetscObject)A,&comm);
2804:   MPI_Comm_size(comm,&size);
2805:   MatDenseGetArrayRead(a->A,&av);
2806:   MatDenseGetArrayRead(b->A,&bv);
2807:   MatDenseGetArrayWrite(c->A,&cv);
2808:   MatDenseGetLDA(a->A,&i);
2809:   PetscBLASIntCast(i,&alda);
2810:   MatDenseGetLDA(b->A,&blda);
2811:   MatDenseGetLDA(c->A,&i);
2812:   PetscBLASIntCast(i,&clda);
2813:   /* copy transpose of B into buf[0] */
2814:   bn      = B->rmap->n;
2815:   sendbuf = abt->buf[0];
2816:   recvbuf = abt->buf[1];
2817:   for (k = 0, j = 0; j < bn; j++) {
2818:     for (i = 0; i < cK; i++, k++) {
2819:       sendbuf[k] = bv[i * blda + j];
2820:     }
2821:   }
2822:   MatDenseRestoreArrayRead(b->A,&bv);
2823:   MPI_Allgatherv(sendbuf, bn * cK, MPIU_SCALAR, recvbuf, abt->recvcounts, abt->recvdispls, MPIU_SCALAR, comm);
2824:   PetscBLASIntCast(cK,&ck);
2825:   PetscBLASIntCast(c->A->rmap->n,&cm);
2826:   PetscBLASIntCast(c->A->cmap->n,&cn);
2827:   if (cm && cn && ck) PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&cm,&cn,&ck,&_DOne,av,&alda,recvbuf,&ck,&_DZero,cv,&clda));
2828:   MatDenseRestoreArrayRead(a->A,&av);
2829:   MatDenseRestoreArrayRead(b->A,&bv);
2830:   MatDenseRestoreArrayWrite(c->A,&cv);
2831:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2832:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2833:   return(0);
2834: }

2836: static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense(Mat A, Mat B, Mat C)
2837: {
2838:   Mat_MatTransMultDense *abt;
2839:   PetscErrorCode        ierr;

2842:   MatCheckProduct(C,3);
2843:   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
2844:   abt = (Mat_MatTransMultDense*)C->product->data;
2845:   switch (abt->alg) {
2846:   case 1:
2847:     MatMatTransposeMultNumeric_MPIDense_MPIDense_Cyclic(A, B, C);
2848:     break;
2849:   default:
2850:     MatMatTransposeMultNumeric_MPIDense_MPIDense_Allgatherv(A, B, C);
2851:     break;
2852:   }
2853:   return(0);
2854: }

2856: PetscErrorCode MatDestroy_MatMatMult_MPIDense_MPIDense(void *data)
2857: {
2858:   PetscErrorCode   ierr;
2859:   Mat_MatMultDense *ab = (Mat_MatMultDense*)data;

2862:   MatDestroy(&ab->Ce);
2863:   MatDestroy(&ab->Ae);
2864:   MatDestroy(&ab->Be);
2865:   PetscFree(ab);
2866:   return(0);
2867: }

2869: #if defined(PETSC_HAVE_ELEMENTAL)
2870: PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C)
2871: {
2872:   PetscErrorCode   ierr;
2873:   Mat_MatMultDense *ab;

2876:   MatCheckProduct(C,3);
2877:   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing product data");
2878:   ab   = (Mat_MatMultDense*)C->product->data;
2879:   MatConvert_MPIDense_Elemental(A,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Ae);
2880:   MatConvert_MPIDense_Elemental(B,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Be);
2881:   MatMatMultNumeric_Elemental(ab->Ae,ab->Be,ab->Ce);
2882:   MatConvert(ab->Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C);
2883:   return(0);
2884: }

2886: static PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
2887: {
2888:   PetscErrorCode   ierr;
2889:   Mat              Ae,Be,Ce;
2890:   Mat_MatMultDense *ab;

2893:   MatCheckProduct(C,4);
2894:   if (C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
2895:   /* check local size of A and B */
2896:   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) {
2897:     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);
2898:   }

2900:   /* create elemental matrices Ae and Be */
2901:   MatCreate(PetscObjectComm((PetscObject)A), &Ae);
2902:   MatSetSizes(Ae,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);
2903:   MatSetType(Ae,MATELEMENTAL);
2904:   MatSetUp(Ae);
2905:   MatSetOption(Ae,MAT_ROW_ORIENTED,PETSC_FALSE);

2907:   MatCreate(PetscObjectComm((PetscObject)B), &Be);
2908:   MatSetSizes(Be,PETSC_DECIDE,PETSC_DECIDE,B->rmap->N,B->cmap->N);
2909:   MatSetType(Be,MATELEMENTAL);
2910:   MatSetUp(Be);
2911:   MatSetOption(Be,MAT_ROW_ORIENTED,PETSC_FALSE);

2913:   /* compute symbolic Ce = Ae*Be */
2914:   MatCreate(PetscObjectComm((PetscObject)C),&Ce);
2915:   MatMatMultSymbolic_Elemental(Ae,Be,fill,Ce);

2917:   /* setup C */
2918:   MatSetSizes(C,A->rmap->n,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
2919:   MatSetType(C,MATDENSE);
2920:   MatSetUp(C);

2922:   /* create data structure for reuse Cdense */
2923:   PetscNew(&ab);
2924:   ab->Ae = Ae;
2925:   ab->Be = Be;
2926:   ab->Ce = Ce;

2928:   C->product->data    = ab;
2929:   C->product->destroy = MatDestroy_MatMatMult_MPIDense_MPIDense;
2930:   C->ops->matmultnumeric = MatMatMultNumeric_MPIDense_MPIDense;
2931:   return(0);
2932: }
2933: #endif
2934: /* ----------------------------------------------- */
2935: #if defined(PETSC_HAVE_ELEMENTAL)
2936: static PetscErrorCode MatProductSetFromOptions_MPIDense_AB(Mat C)
2937: {
2939:   C->ops->matmultsymbolic = MatMatMultSymbolic_MPIDense_MPIDense;
2940:   C->ops->productsymbolic = MatProductSymbolic_AB;
2941:   return(0);
2942: }
2943: #endif

2945: static PetscErrorCode MatProductSetFromOptions_MPIDense_AtB(Mat C)
2946: {
2947:   Mat_Product *product = C->product;
2948:   Mat         A = product->A,B=product->B;

2951:   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend)
2952:     SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
2953:   C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_MPIDense_MPIDense;
2954:   C->ops->productsymbolic = MatProductSymbolic_AtB;
2955:   return(0);
2956: }

2958: static PetscErrorCode MatProductSetFromOptions_MPIDense_ABt(Mat C)
2959: {
2961:   Mat_Product    *product = C->product;
2962:   const char     *algTypes[2] = {"allgatherv","cyclic"};
2963:   PetscInt       alg,nalg = 2;
2964:   PetscBool      flg = PETSC_FALSE;

2967:   /* Set default algorithm */
2968:   alg = 0; /* default is allgatherv */
2969:   PetscStrcmp(product->alg,"default",&flg);
2970:   if (flg) {
2971:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2972:   }

2974:   /* Get runtime option */
2975:   if (product->api_user) {
2976:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatTransposeMult","Mat");
2977:     PetscOptionsEList("-matmattransmult_mpidense_mpidense_via","Algorithmic approach","MatMatTransposeMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2978:     PetscOptionsEnd();
2979:   } else {
2980:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABt","Mat");
2981:     PetscOptionsEList("-matproduct_abt_mpidense_mpidense_via","Algorithmic approach","MatProduct_ABt",algTypes,nalg,algTypes[alg],&alg,&flg);
2982:     PetscOptionsEnd();
2983:   }
2984:   if (flg) {
2985:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2986:   }

2988:   C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_MPIDense_MPIDense;
2989:   C->ops->productsymbolic          = MatProductSymbolic_ABt;
2990:   return(0);
2991: }

2993: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense(Mat C)
2994: {
2996:   Mat_Product    *product = C->product;

2999:   switch (product->type) {
3000: #if defined(PETSC_HAVE_ELEMENTAL)
3001:   case MATPRODUCT_AB:
3002:     MatProductSetFromOptions_MPIDense_AB(C);
3003:     break;
3004: #endif
3005:   case MATPRODUCT_AtB:
3006:     MatProductSetFromOptions_MPIDense_AtB(C);
3007:     break;
3008:   case MATPRODUCT_ABt:
3009:     MatProductSetFromOptions_MPIDense_ABt(C);
3010:     break;
3011:   default:
3012:     break;
3013:   }
3014:   return(0);
3015: }