xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision cc2e6a90c05b27ffec69cb207fe793d447f14420)
1be1d678aSKris Buschelman 
2ed3cc1f0SBarry Smith /*
3ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
4ed3cc1f0SBarry Smith */
5ed3cc1f0SBarry Smith 
604fea9ffSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>    /*I   "petscmat.h"  I*/
88949adfdSHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>
9baa3c1c6SHong Zhang #include <petscblaslapack.h>
108965ea79SLois Curfman McInnes 
11ba8c8a56SBarry Smith #undef __FUNCT__
12ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
13ab92ecdeSBarry Smith /*@
14ab92ecdeSBarry Smith 
15ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
16ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
17ab92ecdeSBarry Smith 
18ab92ecdeSBarry Smith     Input Parameter:
19ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
20ab92ecdeSBarry Smith 
21ab92ecdeSBarry Smith     Output Parameter:
22ab92ecdeSBarry Smith .      B - the inner matrix
23ab92ecdeSBarry Smith 
248e6c10adSSatish Balay     Level: intermediate
258e6c10adSSatish Balay 
26ab92ecdeSBarry Smith @*/
27ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
28ab92ecdeSBarry Smith {
29ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
30ab92ecdeSBarry Smith   PetscErrorCode ierr;
31ace3abfcSBarry Smith   PetscBool      flg;
32ab92ecdeSBarry Smith 
33ab92ecdeSBarry Smith   PetscFunctionBegin;
34251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
352205254eSKarl Rupp   if (flg) *B = mat->A;
362205254eSKarl Rupp   else *B = A;
37ab92ecdeSBarry Smith   PetscFunctionReturn(0);
38ab92ecdeSBarry Smith }
39ab92ecdeSBarry Smith 
40ab92ecdeSBarry Smith #undef __FUNCT__
41ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
42ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
43ba8c8a56SBarry Smith {
44ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
45ba8c8a56SBarry Smith   PetscErrorCode ierr;
46d0f46423SBarry Smith   PetscInt       lrow,rstart = A->rmap->rstart,rend = A->rmap->rend;
47ba8c8a56SBarry Smith 
48ba8c8a56SBarry Smith   PetscFunctionBegin;
49e7e72b3dSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows");
50ba8c8a56SBarry Smith   lrow = row - rstart;
51ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt**)idx,(const PetscScalar**)v);CHKERRQ(ierr);
52ba8c8a56SBarry Smith   PetscFunctionReturn(0);
53ba8c8a56SBarry Smith }
54ba8c8a56SBarry Smith 
55ba8c8a56SBarry Smith #undef __FUNCT__
56ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
57ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
58ba8c8a56SBarry Smith {
59ba8c8a56SBarry Smith   PetscErrorCode ierr;
60ba8c8a56SBarry Smith 
61ba8c8a56SBarry Smith   PetscFunctionBegin;
62ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
63ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
64ba8c8a56SBarry Smith   PetscFunctionReturn(0);
65ba8c8a56SBarry Smith }
66ba8c8a56SBarry Smith 
674a2ae208SSatish Balay #undef __FUNCT__
684a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
6911bd1e4dSLisandro Dalcin PetscErrorCode  MatGetDiagonalBlock_MPIDense(Mat A,Mat *a)
700de54da6SSatish Balay {
710de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
726849ba73SBarry Smith   PetscErrorCode ierr;
73d0f46423SBarry Smith   PetscInt       m = A->rmap->n,rstart = A->rmap->rstart;
7487828ca2SBarry Smith   PetscScalar    *array;
750de54da6SSatish Balay   MPI_Comm       comm;
7611bd1e4dSLisandro Dalcin   Mat            B;
770de54da6SSatish Balay 
780de54da6SSatish Balay   PetscFunctionBegin;
79e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only square matrices supported.");
800de54da6SSatish Balay 
8111bd1e4dSLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);CHKERRQ(ierr);
8211bd1e4dSLisandro Dalcin   if (!B) {
830de54da6SSatish Balay     ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
8411bd1e4dSLisandro Dalcin     ierr = MatCreate(comm,&B);CHKERRQ(ierr);
8511bd1e4dSLisandro Dalcin     ierr = MatSetSizes(B,m,m,m,m);CHKERRQ(ierr);
8611bd1e4dSLisandro Dalcin     ierr = MatSetType(B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr);
878c778c55SBarry Smith     ierr = MatDenseGetArray(mdn->A,&array);CHKERRQ(ierr);
8811bd1e4dSLisandro Dalcin     ierr = MatSeqDenseSetPreallocation(B,array+m*rstart);CHKERRQ(ierr);
898c778c55SBarry Smith     ierr = MatDenseRestoreArray(mdn->A,&array);CHKERRQ(ierr);
9011bd1e4dSLisandro Dalcin     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9111bd1e4dSLisandro Dalcin     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9211bd1e4dSLisandro Dalcin     ierr = PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);CHKERRQ(ierr);
9311bd1e4dSLisandro Dalcin     *a   = B;
9411bd1e4dSLisandro Dalcin     ierr = MatDestroy(&B);CHKERRQ(ierr);
952205254eSKarl Rupp   } else *a = B;
960de54da6SSatish Balay   PetscFunctionReturn(0);
970de54da6SSatish Balay }
980de54da6SSatish Balay 
994a2ae208SSatish Balay #undef __FUNCT__
1004a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
10113f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1028965ea79SLois Curfman McInnes {
10339b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
104dfbe8321SBarry Smith   PetscErrorCode ierr;
105d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
106ace3abfcSBarry Smith   PetscBool      roworiented = A->roworiented;
1078965ea79SLois Curfman McInnes 
1083a40ed3dSBarry Smith   PetscFunctionBegin;
1098965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1105ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
111e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1128965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1138965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
11439b7565bSBarry Smith       if (roworiented) {
11539b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1163a40ed3dSBarry Smith       } else {
1178965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1185ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
119e32f2f54SBarry Smith           if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
12039b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
12139b7565bSBarry Smith         }
1228965ea79SLois Curfman McInnes       }
1232205254eSKarl Rupp     } else if (!A->donotstash) {
1245080c13bSMatthew G Knepley       mat->assembled = PETSC_FALSE;
12539b7565bSBarry Smith       if (roworiented) {
126b400d20cSBarry Smith         ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n,PETSC_FALSE);CHKERRQ(ierr);
127d36fbae8SSatish Balay       } else {
128b400d20cSBarry Smith         ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m,PETSC_FALSE);CHKERRQ(ierr);
12939b7565bSBarry Smith       }
130b49de8d1SLois Curfman McInnes     }
131b49de8d1SLois Curfman McInnes   }
1323a40ed3dSBarry Smith   PetscFunctionReturn(0);
133b49de8d1SLois Curfman McInnes }
134b49de8d1SLois Curfman McInnes 
1354a2ae208SSatish Balay #undef __FUNCT__
1364a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
13713f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
138b49de8d1SLois Curfman McInnes {
139b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
140dfbe8321SBarry Smith   PetscErrorCode ierr;
141d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
142b49de8d1SLois Curfman McInnes 
1433a40ed3dSBarry Smith   PetscFunctionBegin;
144b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
145e32f2f54SBarry Smith     if (idxm[i] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
146e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
147b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
148b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
149b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
150e32f2f54SBarry Smith         if (idxn[j] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
151e7e72b3dSBarry Smith         if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
152b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
153b49de8d1SLois Curfman McInnes       }
154e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
1558965ea79SLois Curfman McInnes   }
1563a40ed3dSBarry Smith   PetscFunctionReturn(0);
1578965ea79SLois Curfman McInnes }
1588965ea79SLois Curfman McInnes 
1594a2ae208SSatish Balay #undef __FUNCT__
1608c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_MPIDense"
1618c778c55SBarry Smith PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[])
162ff14e315SSatish Balay {
163ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
164dfbe8321SBarry Smith   PetscErrorCode ierr;
165ff14e315SSatish Balay 
1663a40ed3dSBarry Smith   PetscFunctionBegin;
1678c778c55SBarry Smith   ierr = MatDenseGetArray(a->A,array);CHKERRQ(ierr);
1683a40ed3dSBarry Smith   PetscFunctionReturn(0);
169ff14e315SSatish Balay }
170ff14e315SSatish Balay 
1714a2ae208SSatish Balay #undef __FUNCT__
1724a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
1734aa3045dSJed Brown static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
174ca3fa75bSLois Curfman McInnes {
175ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat  = (Mat_MPIDense*)A->data,*newmatd;
176ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
1776849ba73SBarry Smith   PetscErrorCode ierr;
1784aa3045dSJed Brown   PetscInt       i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols;
1795d0c19d7SBarry Smith   const PetscInt *irow,*icol;
18087828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
181ca3fa75bSLois Curfman McInnes   Mat            newmat;
1824aa3045dSJed Brown   IS             iscol_local;
183ca3fa75bSLois Curfman McInnes 
184ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
1854aa3045dSJed Brown   ierr = ISAllGather(iscol,&iscol_local);CHKERRQ(ierr);
186ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
1874aa3045dSJed Brown   ierr = ISGetIndices(iscol_local,&icol);CHKERRQ(ierr);
188b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
189b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
1904aa3045dSJed Brown   ierr = ISGetSize(iscol,&Ncols);CHKERRQ(ierr); /* global number of columns, size of iscol_local */
191ca3fa75bSLois Curfman McInnes 
192ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1937eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1947eba5e9cSLois Curfman McInnes      original matrix! */
195ca3fa75bSLois Curfman McInnes 
196ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1977eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
198ca3fa75bSLois Curfman McInnes 
199ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
200ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
201e32f2f54SBarry Smith     /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2027eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
203ca3fa75bSLois Curfman McInnes     newmat = *B;
204ca3fa75bSLois Curfman McInnes   } else {
205ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
206ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr);
2074aa3045dSJed Brown     ierr = MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);CHKERRQ(ierr);
2087adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
2090298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
210ca3fa75bSLois Curfman McInnes   }
211ca3fa75bSLois Curfman McInnes 
212ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
213ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
214ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense*)newmatd->A->data)->v;
215ca3fa75bSLois Curfman McInnes 
2164aa3045dSJed Brown   for (i=0; i<Ncols; i++) {
21725a33276SHong Zhang     av = v + ((Mat_SeqDense*)mat->A->data)->lda*icol[i];
218ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2197eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
220ca3fa75bSLois Curfman McInnes     }
221ca3fa75bSLois Curfman McInnes   }
222ca3fa75bSLois Curfman McInnes 
223ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
224ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
225ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
226ca3fa75bSLois Curfman McInnes 
227ca3fa75bSLois Curfman McInnes   /* Free work space */
228ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
2295bdf786aSShri Abhyankar   ierr = ISRestoreIndices(iscol_local,&icol);CHKERRQ(ierr);
23032bb1f2dSLisandro Dalcin   ierr = ISDestroy(&iscol_local);CHKERRQ(ierr);
231ca3fa75bSLois Curfman McInnes   *B   = newmat;
232ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
233ca3fa75bSLois Curfman McInnes }
234ca3fa75bSLois Curfman McInnes 
2354a2ae208SSatish Balay #undef __FUNCT__
2368c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_MPIDense"
2378c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[])
238ff14e315SSatish Balay {
23973a71a0fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
24073a71a0fSBarry Smith   PetscErrorCode ierr;
24173a71a0fSBarry Smith 
2423a40ed3dSBarry Smith   PetscFunctionBegin;
2438c778c55SBarry Smith   ierr = MatDenseRestoreArray(a->A,array);CHKERRQ(ierr);
2443a40ed3dSBarry Smith   PetscFunctionReturn(0);
245ff14e315SSatish Balay }
246ff14e315SSatish Balay 
2474a2ae208SSatish Balay #undef __FUNCT__
2484a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
249dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2508965ea79SLois Curfman McInnes {
25139ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
252ce94432eSBarry Smith   MPI_Comm       comm;
253dfbe8321SBarry Smith   PetscErrorCode ierr;
25413f74950SBarry Smith   PetscInt       nstash,reallocs;
2558965ea79SLois Curfman McInnes   InsertMode     addv;
2568965ea79SLois Curfman McInnes 
2573a40ed3dSBarry Smith   PetscFunctionBegin;
258ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
2598965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
260b2566f29SBarry Smith   ierr = MPIU_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);CHKERRQ(ierr);
261ce94432eSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
262e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2638965ea79SLois Curfman McInnes 
264d0f46423SBarry Smith   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr);
2658798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
266ae15b995SBarry Smith   ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
2673a40ed3dSBarry Smith   PetscFunctionReturn(0);
2688965ea79SLois Curfman McInnes }
2698965ea79SLois Curfman McInnes 
2704a2ae208SSatish Balay #undef __FUNCT__
2714a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
272dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2738965ea79SLois Curfman McInnes {
27439ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn=(Mat_MPIDense*)mat->data;
2756849ba73SBarry Smith   PetscErrorCode ierr;
27613f74950SBarry Smith   PetscInt       i,*row,*col,flg,j,rstart,ncols;
27713f74950SBarry Smith   PetscMPIInt    n;
27887828ca2SBarry Smith   PetscScalar    *val;
2798965ea79SLois Curfman McInnes 
2803a40ed3dSBarry Smith   PetscFunctionBegin;
2818965ea79SLois Curfman McInnes   /*  wait on receives */
2827ef1d9bdSSatish Balay   while (1) {
2838798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2847ef1d9bdSSatish Balay     if (!flg) break;
2858965ea79SLois Curfman McInnes 
2867ef1d9bdSSatish Balay     for (i=0; i<n;) {
2877ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2882205254eSKarl Rupp       for (j=i,rstart=row[j]; j<n; j++) {
2892205254eSKarl Rupp         if (row[j] != rstart) break;
2902205254eSKarl Rupp       }
2917ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2927ef1d9bdSSatish Balay       else       ncols = n-i;
2937ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2944b4eb8d3SJed Brown       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);CHKERRQ(ierr);
2957ef1d9bdSSatish Balay       i    = j;
2968965ea79SLois Curfman McInnes     }
2977ef1d9bdSSatish Balay   }
2988798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2998965ea79SLois Curfman McInnes 
30039ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
30139ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3028965ea79SLois Curfman McInnes 
3036d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
30439ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3058965ea79SLois Curfman McInnes   }
3063a40ed3dSBarry Smith   PetscFunctionReturn(0);
3078965ea79SLois Curfman McInnes }
3088965ea79SLois Curfman McInnes 
3094a2ae208SSatish Balay #undef __FUNCT__
3104a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
311dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3128965ea79SLois Curfman McInnes {
313dfbe8321SBarry Smith   PetscErrorCode ierr;
31439ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3153a40ed3dSBarry Smith 
3163a40ed3dSBarry Smith   PetscFunctionBegin;
3173a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3183a40ed3dSBarry Smith   PetscFunctionReturn(0);
3198965ea79SLois Curfman McInnes }
3208965ea79SLois Curfman McInnes 
3218965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3228965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3238965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3243501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3258965ea79SLois Curfman McInnes    routine.
3268965ea79SLois Curfman McInnes */
3274a2ae208SSatish Balay #undef __FUNCT__
3284a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
3292b40b63fSBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
3308965ea79SLois Curfman McInnes {
33139ddd567SLois Curfman McInnes   Mat_MPIDense      *l = (Mat_MPIDense*)A->data;
3326849ba73SBarry Smith   PetscErrorCode    ierr;
333d0f46423SBarry Smith   PetscInt          i,*owners = A->rmap->range;
33476ec1555SBarry Smith   PetscInt          *sizes,j,idx,nsends;
33513f74950SBarry Smith   PetscInt          nmax,*svalues,*starts,*owner,nrecvs;
3367adad957SLisandro Dalcin   PetscInt          *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source;
33713f74950SBarry Smith   PetscInt          *lens,*lrows,*values;
33813f74950SBarry Smith   PetscMPIInt       n,imdex,rank = l->rank,size = l->size;
339ce94432eSBarry Smith   MPI_Comm          comm;
3408965ea79SLois Curfman McInnes   MPI_Request       *send_waits,*recv_waits;
3418965ea79SLois Curfman McInnes   MPI_Status        recv_status,*send_status;
342ace3abfcSBarry Smith   PetscBool         found;
34397b48c8fSBarry Smith   const PetscScalar *xx;
34497b48c8fSBarry Smith   PetscScalar       *bb;
3458965ea79SLois Curfman McInnes 
3463a40ed3dSBarry Smith   PetscFunctionBegin;
347ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
348ce94432eSBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(comm,PETSC_ERR_SUP,"Only handles square matrices");
349b9679d65SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership");
3508965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
351f628708eSJed Brown   ierr = PetscCalloc1(2*size,&sizes);CHKERRQ(ierr);
352037dbc42SBarry Smith   ierr = PetscMalloc1(N+1,&owner);CHKERRQ(ierr);  /* see note*/
3538965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3548965ea79SLois Curfman McInnes     idx   = rows[i];
35535d8aa7fSBarry Smith     found = PETSC_FALSE;
3568965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3578965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
35876ec1555SBarry Smith         sizes[2*j]++; sizes[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3598965ea79SLois Curfman McInnes       }
3608965ea79SLois Curfman McInnes     }
361e32f2f54SBarry Smith     if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3628965ea79SLois Curfman McInnes   }
3632205254eSKarl Rupp   nsends = 0;
36476ec1555SBarry Smith   for (i=0; i<size; i++) nsends += sizes[2*i+1];
3658965ea79SLois Curfman McInnes 
3668965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
36776ec1555SBarry Smith   ierr = PetscMaxSum(comm,sizes,&nmax,&nrecvs);CHKERRQ(ierr);
3688965ea79SLois Curfman McInnes 
3698965ea79SLois Curfman McInnes   /* post receives:   */
370785e854fSJed Brown   ierr = PetscMalloc1((nrecvs+1)*(nmax+1),&rvalues);CHKERRQ(ierr);
371854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
3728965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
37313f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3748965ea79SLois Curfman McInnes   }
3758965ea79SLois Curfman McInnes 
3768965ea79SLois Curfman McInnes   /* do sends:
3778965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3788965ea79SLois Curfman McInnes          the ith processor
3798965ea79SLois Curfman McInnes   */
380854ce69bSBarry Smith   ierr = PetscMalloc1(N+1,&svalues);CHKERRQ(ierr);
381854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
382854ce69bSBarry Smith   ierr = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
3838965ea79SLois Curfman McInnes 
3848965ea79SLois Curfman McInnes   starts[0] = 0;
38576ec1555SBarry Smith   for (i=1; i<size; i++) starts[i] = starts[i-1] + sizes[2*i-2];
3862205254eSKarl Rupp   for (i=0; i<N; i++) svalues[starts[owner[i]]++] = rows[i];
3872205254eSKarl Rupp 
3882205254eSKarl Rupp   starts[0] = 0;
38976ec1555SBarry Smith   for (i=1; i<size+1; i++) starts[i] = starts[i-1] + sizes[2*i-2];
3908965ea79SLois Curfman McInnes   count = 0;
3918965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
39276ec1555SBarry Smith     if (sizes[2*i+1]) {
39376ec1555SBarry Smith       ierr = MPI_Isend(svalues+starts[i],sizes[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3948965ea79SLois Curfman McInnes     }
3958965ea79SLois Curfman McInnes   }
396606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3978965ea79SLois Curfman McInnes 
3988965ea79SLois Curfman McInnes   base = owners[rank];
3998965ea79SLois Curfman McInnes 
4008965ea79SLois Curfman McInnes   /*  wait on receives */
401dcca6d9dSJed Brown   ierr  = PetscMalloc2(nrecvs,&lens,nrecvs,&source);CHKERRQ(ierr);
40274ed9c26SBarry Smith   count = nrecvs;
40374ed9c26SBarry Smith   slen  = 0;
4048965ea79SLois Curfman McInnes   while (count) {
405ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4068965ea79SLois Curfman McInnes     /* unpack receives into our local space */
40713f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4082205254eSKarl Rupp 
4098965ea79SLois Curfman McInnes     source[imdex] = recv_status.MPI_SOURCE;
4108965ea79SLois Curfman McInnes     lens[imdex]   = n;
4118965ea79SLois Curfman McInnes     slen += n;
4128965ea79SLois Curfman McInnes     count--;
4138965ea79SLois Curfman McInnes   }
414606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4158965ea79SLois Curfman McInnes 
4168965ea79SLois Curfman McInnes   /* move the data into the send scatter */
417854ce69bSBarry Smith   ierr  = PetscMalloc1(slen+1,&lrows);CHKERRQ(ierr);
4188965ea79SLois Curfman McInnes   count = 0;
4198965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4208965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4218965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4228965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4238965ea79SLois Curfman McInnes     }
4248965ea79SLois Curfman McInnes   }
425606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
42674ed9c26SBarry Smith   ierr = PetscFree2(lens,source);CHKERRQ(ierr);
427606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
42876ec1555SBarry Smith   ierr = PetscFree(sizes);CHKERRQ(ierr);
4298965ea79SLois Curfman McInnes 
43097b48c8fSBarry Smith   /* fix right hand side if needed */
43197b48c8fSBarry Smith   if (x && b) {
43297b48c8fSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
43397b48c8fSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
43497b48c8fSBarry Smith     for (i=0; i<slen; i++) {
43597b48c8fSBarry Smith       bb[lrows[i]] = diag*xx[lrows[i]];
43697b48c8fSBarry Smith     }
43797b48c8fSBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
43897b48c8fSBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
43997b48c8fSBarry Smith   }
44097b48c8fSBarry Smith 
4418965ea79SLois Curfman McInnes   /* actually zap the local rows */
442b9679d65SBarry Smith   ierr = MatZeroRows(l->A,slen,lrows,0.0,0,0);CHKERRQ(ierr);
443e2eb51b1SBarry Smith   if (diag != 0.0) {
444b9679d65SBarry Smith     Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data;
445b9679d65SBarry Smith     PetscInt     m   = ll->lda, i;
446b9679d65SBarry Smith 
447b9679d65SBarry Smith     for (i=0; i<slen; i++) {
448b9679d65SBarry Smith       ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag;
449b9679d65SBarry Smith     }
450b9679d65SBarry Smith   }
451606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4528965ea79SLois Curfman McInnes 
4538965ea79SLois Curfman McInnes   /* wait on sends */
4548965ea79SLois Curfman McInnes   if (nsends) {
455785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
456ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
457606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4588965ea79SLois Curfman McInnes   }
459606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
460606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4613a40ed3dSBarry Smith   PetscFunctionReturn(0);
4628965ea79SLois Curfman McInnes }
4638965ea79SLois Curfman McInnes 
464*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
465*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
466*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
467*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);
468*cc2e6a90SBarry Smith 
4694a2ae208SSatish Balay #undef __FUNCT__
4704a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
471dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4728965ea79SLois Curfman McInnes {
47339ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
474dfbe8321SBarry Smith   PetscErrorCode ierr;
475c456f294SBarry Smith 
4763a40ed3dSBarry Smith   PetscFunctionBegin;
477ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
478ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
479*cc2e6a90SBarry Smith   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4803a40ed3dSBarry Smith   PetscFunctionReturn(0);
4818965ea79SLois Curfman McInnes }
4828965ea79SLois Curfman McInnes 
4834a2ae208SSatish Balay #undef __FUNCT__
4844a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
485dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4868965ea79SLois Curfman McInnes {
48739ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
488dfbe8321SBarry Smith   PetscErrorCode ierr;
489c456f294SBarry Smith 
4903a40ed3dSBarry Smith   PetscFunctionBegin;
491ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
492ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
493*cc2e6a90SBarry Smith   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4943a40ed3dSBarry Smith   PetscFunctionReturn(0);
4958965ea79SLois Curfman McInnes }
4968965ea79SLois Curfman McInnes 
4974a2ae208SSatish Balay #undef __FUNCT__
4984a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
499dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
500096963f5SLois Curfman McInnes {
501096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
502dfbe8321SBarry Smith   PetscErrorCode ierr;
50387828ca2SBarry Smith   PetscScalar    zero = 0.0;
504096963f5SLois Curfman McInnes 
5053a40ed3dSBarry Smith   PetscFunctionBegin;
5062dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
507*cc2e6a90SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
508ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
509ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5103a40ed3dSBarry Smith   PetscFunctionReturn(0);
511096963f5SLois Curfman McInnes }
512096963f5SLois Curfman McInnes 
5134a2ae208SSatish Balay #undef __FUNCT__
5144a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
515dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
516096963f5SLois Curfman McInnes {
517096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
518dfbe8321SBarry Smith   PetscErrorCode ierr;
519096963f5SLois Curfman McInnes 
5203a40ed3dSBarry Smith   PetscFunctionBegin;
5213501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
522*cc2e6a90SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
523ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
524ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5253a40ed3dSBarry Smith   PetscFunctionReturn(0);
526096963f5SLois Curfman McInnes }
527096963f5SLois Curfman McInnes 
5284a2ae208SSatish Balay #undef __FUNCT__
5294a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
530dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5318965ea79SLois Curfman McInnes {
53239ddd567SLois Curfman McInnes   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
533096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
534dfbe8321SBarry Smith   PetscErrorCode ierr;
535d0f46423SBarry Smith   PetscInt       len,i,n,m = A->rmap->n,radd;
53687828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
537ed3cc1f0SBarry Smith 
5383a40ed3dSBarry Smith   PetscFunctionBegin;
5392dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5401ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
541096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
542e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
543d0f46423SBarry Smith   len  = PetscMin(a->A->rmap->n,a->A->cmap->n);
544d0f46423SBarry Smith   radd = A->rmap->rstart*m;
54544cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
546096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
547096963f5SLois Curfman McInnes   }
5481ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5493a40ed3dSBarry Smith   PetscFunctionReturn(0);
5508965ea79SLois Curfman McInnes }
5518965ea79SLois Curfman McInnes 
5524a2ae208SSatish Balay #undef __FUNCT__
5534a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
554dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5558965ea79SLois Curfman McInnes {
5563501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
557dfbe8321SBarry Smith   PetscErrorCode ierr;
558ed3cc1f0SBarry Smith 
5593a40ed3dSBarry Smith   PetscFunctionBegin;
560aa482453SBarry Smith #if defined(PETSC_USE_LOG)
561d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N);
5628965ea79SLois Curfman McInnes #endif
5638798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
5646bf464f9SBarry Smith   ierr = MatDestroy(&mdn->A);CHKERRQ(ierr);
5656bf464f9SBarry Smith   ierr = VecDestroy(&mdn->lvec);CHKERRQ(ierr);
5666bf464f9SBarry Smith   ierr = VecScatterDestroy(&mdn->Mvctx);CHKERRQ(ierr);
56701b82886SBarry Smith 
568bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
569dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
5708baccfbdSHong Zhang 
5718baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);CHKERRQ(ierr);
5728baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);CHKERRQ(ierr);
5738baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
5748baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",NULL);CHKERRQ(ierr);
5758baccfbdSHong Zhang #endif
576bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",NULL);CHKERRQ(ierr);
577bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",NULL);CHKERRQ(ierr);
578bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
579bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
580bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5818baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5828baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5838baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5843a40ed3dSBarry Smith   PetscFunctionReturn(0);
5858965ea79SLois Curfman McInnes }
58639ddd567SLois Curfman McInnes 
5874a2ae208SSatish Balay #undef __FUNCT__
5884a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5896849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5908965ea79SLois Curfman McInnes {
59139ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
592dfbe8321SBarry Smith   PetscErrorCode    ierr;
593aa05aa95SBarry Smith   PetscViewerFormat format;
594aa05aa95SBarry Smith   int               fd;
595d0f46423SBarry Smith   PetscInt          header[4],mmax,N = mat->cmap->N,i,j,m,k;
596aa05aa95SBarry Smith   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
597578230a0SSatish Balay   PetscScalar       *work,*v,*vv;
598aa05aa95SBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;
5997056b6fcSBarry Smith 
6003a40ed3dSBarry Smith   PetscFunctionBegin;
60139ddd567SLois Curfman McInnes   if (mdn->size == 1) {
60239ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
603aa05aa95SBarry Smith   } else {
604aa05aa95SBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
605ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr);
606ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
607aa05aa95SBarry Smith 
608aa05aa95SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
609f4403165SShri Abhyankar     if (format == PETSC_VIEWER_NATIVE) {
610aa05aa95SBarry Smith 
611aa05aa95SBarry Smith       if (!rank) {
612aa05aa95SBarry Smith         /* store the matrix as a dense matrix */
6130700a824SBarry Smith         header[0] = MAT_FILE_CLASSID;
614d0f46423SBarry Smith         header[1] = mat->rmap->N;
615aa05aa95SBarry Smith         header[2] = N;
616aa05aa95SBarry Smith         header[3] = MATRIX_BINARY_FORMAT_DENSE;
617aa05aa95SBarry Smith         ierr      = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
618aa05aa95SBarry Smith 
619aa05aa95SBarry Smith         /* get largest work array needed for transposing array */
620d0f46423SBarry Smith         mmax = mat->rmap->n;
621aa05aa95SBarry Smith         for (i=1; i<size; i++) {
622d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
6238965ea79SLois Curfman McInnes         }
624785e854fSJed Brown         ierr = PetscMalloc1(mmax*N,&work);CHKERRQ(ierr);
625aa05aa95SBarry Smith 
626aa05aa95SBarry Smith         /* write out local array, by rows */
627d0f46423SBarry Smith         m = mat->rmap->n;
628aa05aa95SBarry Smith         v = a->v;
629aa05aa95SBarry Smith         for (j=0; j<N; j++) {
630aa05aa95SBarry Smith           for (i=0; i<m; i++) {
631578230a0SSatish Balay             work[j + i*N] = *v++;
632aa05aa95SBarry Smith           }
633aa05aa95SBarry Smith         }
634aa05aa95SBarry Smith         ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
635aa05aa95SBarry Smith         /* get largest work array to receive messages from other processes, excludes process zero */
636aa05aa95SBarry Smith         mmax = 0;
637aa05aa95SBarry Smith         for (i=1; i<size; i++) {
638d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
639aa05aa95SBarry Smith         }
640785e854fSJed Brown         ierr = PetscMalloc1(mmax*N,&vv);CHKERRQ(ierr);
641aa05aa95SBarry Smith         for (k = 1; k < size; k++) {
642f8009846SMatthew Knepley           v    = vv;
643d0f46423SBarry Smith           m    = mat->rmap->range[k+1] - mat->rmap->range[k];
644ce94432eSBarry Smith           ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
645aa05aa95SBarry Smith 
646aa05aa95SBarry Smith           for (j = 0; j < N; j++) {
647aa05aa95SBarry Smith             for (i = 0; i < m; i++) {
648578230a0SSatish Balay               work[j + i*N] = *v++;
649aa05aa95SBarry Smith             }
650aa05aa95SBarry Smith           }
651aa05aa95SBarry Smith           ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
652aa05aa95SBarry Smith         }
653aa05aa95SBarry Smith         ierr = PetscFree(work);CHKERRQ(ierr);
654578230a0SSatish Balay         ierr = PetscFree(vv);CHKERRQ(ierr);
655aa05aa95SBarry Smith       } else {
656ce94432eSBarry Smith         ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
657aa05aa95SBarry Smith       }
6586a9046bcSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE)");
659aa05aa95SBarry Smith   }
6603a40ed3dSBarry Smith   PetscFunctionReturn(0);
6618965ea79SLois Curfman McInnes }
6628965ea79SLois Curfman McInnes 
6637da1fb6eSBarry Smith extern PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
6649804daf3SBarry Smith #include <petscdraw.h>
6654a2ae208SSatish Balay #undef __FUNCT__
6664a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
6676849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6688965ea79SLois Curfman McInnes {
66939ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
670dfbe8321SBarry Smith   PetscErrorCode    ierr;
6717da1fb6eSBarry Smith   PetscMPIInt       rank = mdn->rank;
67219fd82e9SBarry Smith   PetscViewerType   vtype;
673ace3abfcSBarry Smith   PetscBool         iascii,isdraw;
674b0a32e0cSBarry Smith   PetscViewer       sviewer;
675f3ef73ceSBarry Smith   PetscViewerFormat format;
6768965ea79SLois Curfman McInnes 
6773a40ed3dSBarry Smith   PetscFunctionBegin;
678251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
679251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
68032077d6dSBarry Smith   if (iascii) {
681b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
682b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
683456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6844e220ebcSLois Curfman McInnes       MatInfo info;
685888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
6861575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
6877b23a99aSBarry Smith       ierr = 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);CHKERRQ(ierr);
688b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6891575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
6905d00a290SHong Zhang       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
6913a40ed3dSBarry Smith       PetscFunctionReturn(0);
692fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6933a40ed3dSBarry Smith       PetscFunctionReturn(0);
6948965ea79SLois Curfman McInnes     }
695f1af5d2fSBarry Smith   } else if (isdraw) {
696b0a32e0cSBarry Smith     PetscDraw draw;
697ace3abfcSBarry Smith     PetscBool isnull;
698f1af5d2fSBarry Smith 
699b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
700b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
701f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
702f1af5d2fSBarry Smith   }
70377ed5343SBarry Smith 
7047da1fb6eSBarry Smith   {
7058965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
7068965ea79SLois Curfman McInnes     Mat         A;
707d0f46423SBarry Smith     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
708ba8c8a56SBarry Smith     PetscInt    *cols;
709ba8c8a56SBarry Smith     PetscScalar *vals;
7108965ea79SLois Curfman McInnes 
711ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr);
7128965ea79SLois Curfman McInnes     if (!rank) {
713f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
7143a40ed3dSBarry Smith     } else {
715f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
7168965ea79SLois Curfman McInnes     }
7177adad957SLisandro Dalcin     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
718878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
7190298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(A,NULL);CHKERRQ(ierr);
7203bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)A);CHKERRQ(ierr);
7218965ea79SLois Curfman McInnes 
72239ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
72339ddd567SLois Curfman McInnes        but it's quick for now */
72451022da4SBarry Smith     A->insertmode = INSERT_VALUES;
7252205254eSKarl Rupp 
7262205254eSKarl Rupp     row = mat->rmap->rstart;
7272205254eSKarl Rupp     m   = mdn->A->rmap->n;
7288965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
729ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
730ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
731ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
73239ddd567SLois Curfman McInnes       row++;
7338965ea79SLois Curfman McInnes     }
7348965ea79SLois Curfman McInnes 
7356d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7366d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7373f08860eSBarry Smith     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
738b9b97703SBarry Smith     if (!rank) {
7397da1fb6eSBarry Smith         ierr = MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
7408965ea79SLois Curfman McInnes     }
7413f08860eSBarry Smith     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
742b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7436bf464f9SBarry Smith     ierr = MatDestroy(&A);CHKERRQ(ierr);
7448965ea79SLois Curfman McInnes   }
7453a40ed3dSBarry Smith   PetscFunctionReturn(0);
7468965ea79SLois Curfman McInnes }
7478965ea79SLois Curfman McInnes 
7484a2ae208SSatish Balay #undef __FUNCT__
7494a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
750dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
7518965ea79SLois Curfman McInnes {
752dfbe8321SBarry Smith   PetscErrorCode ierr;
753ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw,issocket;
7548965ea79SLois Curfman McInnes 
755433994e6SBarry Smith   PetscFunctionBegin;
756251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
757251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
758251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr);
759251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
7600f5bd95cSBarry Smith 
76132077d6dSBarry Smith   if (iascii || issocket || isdraw) {
762f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7630f5bd95cSBarry Smith   } else if (isbinary) {
7643a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
76511aeaf0aSBarry Smith   }
7663a40ed3dSBarry Smith   PetscFunctionReturn(0);
7678965ea79SLois Curfman McInnes }
7688965ea79SLois Curfman McInnes 
7694a2ae208SSatish Balay #undef __FUNCT__
7704a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
771dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7728965ea79SLois Curfman McInnes {
7733501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7743501a2bdSLois Curfman McInnes   Mat            mdn  = mat->A;
775dfbe8321SBarry Smith   PetscErrorCode ierr;
776329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7778965ea79SLois Curfman McInnes 
7783a40ed3dSBarry Smith   PetscFunctionBegin;
7794e220ebcSLois Curfman McInnes   info->block_size = 1.0;
7802205254eSKarl Rupp 
7814e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7822205254eSKarl Rupp 
7834e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7844e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7858965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7864e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7874e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7884e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7894e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7904e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7918965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
792b2566f29SBarry Smith     ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
7932205254eSKarl Rupp 
7944e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7954e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7964e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7974e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7984e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7998965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
800b2566f29SBarry Smith     ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
8012205254eSKarl Rupp 
8024e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
8034e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
8044e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
8054e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
8064e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
8078965ea79SLois Curfman McInnes   }
8084e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
8094e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
8104e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
8113a40ed3dSBarry Smith   PetscFunctionReturn(0);
8128965ea79SLois Curfman McInnes }
8138965ea79SLois Curfman McInnes 
8144a2ae208SSatish Balay #undef __FUNCT__
8154a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
816ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg)
8178965ea79SLois Curfman McInnes {
81839ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
819dfbe8321SBarry Smith   PetscErrorCode ierr;
8208965ea79SLois Curfman McInnes 
8213a40ed3dSBarry Smith   PetscFunctionBegin;
82212c028f9SKris Buschelman   switch (op) {
823512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
82412c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
82512c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
82643674050SBarry Smith     MatCheckPreallocated(A,1);
8274e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82812c028f9SKris Buschelman     break;
82912c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
83043674050SBarry Smith     MatCheckPreallocated(A,1);
8314e0d8c25SBarry Smith     a->roworiented = flg;
8324e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
83312c028f9SKris Buschelman     break;
8344e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
83513fa8e87SLisandro Dalcin   case MAT_KEEP_NONZERO_PATTERN:
83612c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
837290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
83812c028f9SKris Buschelman     break;
83912c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
8404e0d8c25SBarry Smith     a->donotstash = flg;
84112c028f9SKris Buschelman     break;
84277e54ba9SKris Buschelman   case MAT_SYMMETRIC:
84377e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8449a4540c5SBarry Smith   case MAT_HERMITIAN:
8459a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
846600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
847290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
84877e54ba9SKris Buschelman     break;
84912c028f9SKris Buschelman   default:
850e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
8513a40ed3dSBarry Smith   }
8523a40ed3dSBarry Smith   PetscFunctionReturn(0);
8538965ea79SLois Curfman McInnes }
8548965ea79SLois Curfman McInnes 
8558965ea79SLois Curfman McInnes 
8564a2ae208SSatish Balay #undef __FUNCT__
8574a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
858dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8595b2fa520SLois Curfman McInnes {
8605b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8615b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
86287828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
863dfbe8321SBarry Smith   PetscErrorCode ierr;
864d0f46423SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n;
8655b2fa520SLois Curfman McInnes 
8665b2fa520SLois Curfman McInnes   PetscFunctionBegin;
86772d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8685b2fa520SLois Curfman McInnes   if (ll) {
86972d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
870e32f2f54SBarry Smith     if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8711ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8725b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8735b2fa520SLois Curfman McInnes       x = l[i];
8745b2fa520SLois Curfman McInnes       v = mat->v + i;
8755b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8765b2fa520SLois Curfman McInnes     }
8771ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
878efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8795b2fa520SLois Curfman McInnes   }
8805b2fa520SLois Curfman McInnes   if (rr) {
881175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
882e32f2f54SBarry Smith     if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
883ca9f406cSSatish Balay     ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
884ca9f406cSSatish Balay     ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8851ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8865b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8875b2fa520SLois Curfman McInnes       x = r[i];
8885b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8892205254eSKarl Rupp       for (j=0; j<m; j++) (*v++) *= x;
8905b2fa520SLois Curfman McInnes     }
8911ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
892efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8935b2fa520SLois Curfman McInnes   }
8945b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8955b2fa520SLois Curfman McInnes }
8965b2fa520SLois Curfman McInnes 
8974a2ae208SSatish Balay #undef __FUNCT__
8984a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
899dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
900096963f5SLois Curfman McInnes {
9013501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
9023501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
903dfbe8321SBarry Smith   PetscErrorCode ierr;
90413f74950SBarry Smith   PetscInt       i,j;
905329f5518SBarry Smith   PetscReal      sum = 0.0;
90687828ca2SBarry Smith   PetscScalar    *v  = mat->v;
9073501a2bdSLois Curfman McInnes 
9083a40ed3dSBarry Smith   PetscFunctionBegin;
9093501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
910064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
9113501a2bdSLois Curfman McInnes   } else {
9123501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
913d0f46423SBarry Smith       for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) {
914329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
9153501a2bdSLois Curfman McInnes       }
916b2566f29SBarry Smith       ierr = MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
9178f1a2a5eSBarry Smith       *nrm = PetscSqrtReal(*nrm);
918dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr);
9193a40ed3dSBarry Smith     } else if (type == NORM_1) {
920329f5518SBarry Smith       PetscReal *tmp,*tmp2;
921dcca6d9dSJed Brown       ierr = PetscMalloc2(A->cmap->N,&tmp,A->cmap->N,&tmp2);CHKERRQ(ierr);
92274ed9c26SBarry Smith       ierr = PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
92374ed9c26SBarry Smith       ierr = PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
924064f8208SBarry Smith       *nrm = 0.0;
9253501a2bdSLois Curfman McInnes       v    = mat->v;
926d0f46423SBarry Smith       for (j=0; j<mdn->A->cmap->n; j++) {
927d0f46423SBarry Smith         for (i=0; i<mdn->A->rmap->n; i++) {
92867e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
9293501a2bdSLois Curfman McInnes         }
9303501a2bdSLois Curfman McInnes       }
931b2566f29SBarry Smith       ierr = MPIU_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
932d0f46423SBarry Smith       for (j=0; j<A->cmap->N; j++) {
933064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
9343501a2bdSLois Curfman McInnes       }
9358627564fSBarry Smith       ierr = PetscFree2(tmp,tmp2);CHKERRQ(ierr);
936d0f46423SBarry Smith       ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
9373a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
938329f5518SBarry Smith       PetscReal ntemp;
9393501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
940b2566f29SBarry Smith       ierr = MPIU_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
941ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm");
9423501a2bdSLois Curfman McInnes   }
9433a40ed3dSBarry Smith   PetscFunctionReturn(0);
9443501a2bdSLois Curfman McInnes }
9453501a2bdSLois Curfman McInnes 
9464a2ae208SSatish Balay #undef __FUNCT__
9474a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
948fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
9493501a2bdSLois Curfman McInnes {
9503501a2bdSLois Curfman McInnes   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
9513501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
9523501a2bdSLois Curfman McInnes   Mat            B;
953d0f46423SBarry Smith   PetscInt       M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart;
9546849ba73SBarry Smith   PetscErrorCode ierr;
95513f74950SBarry Smith   PetscInt       j,i;
95687828ca2SBarry Smith   PetscScalar    *v;
9573501a2bdSLois Curfman McInnes 
9583a40ed3dSBarry Smith   PetscFunctionBegin;
959ce94432eSBarry Smith   if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place");
960fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
961ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
962d0f46423SBarry Smith     ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr);
9637adad957SLisandro Dalcin     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
9640298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr);
965fc4dec0aSBarry Smith   } else {
966fc4dec0aSBarry Smith     B = *matout;
967fc4dec0aSBarry Smith   }
9683501a2bdSLois Curfman McInnes 
969d0f46423SBarry Smith   m    = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v;
970785e854fSJed Brown   ierr = PetscMalloc1(m,&rwork);CHKERRQ(ierr);
9713501a2bdSLois Curfman McInnes   for (i=0; i<m; i++) rwork[i] = rstart + i;
9721acff37aSSatish Balay   for (j=0; j<n; j++) {
9733501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9743501a2bdSLois Curfman McInnes     v   += m;
9753501a2bdSLois Curfman McInnes   }
976606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9776d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9786d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
979815cbec1SBarry Smith   if (reuse == MAT_INITIAL_MATRIX || *matout != A) {
9803501a2bdSLois Curfman McInnes     *matout = B;
9813501a2bdSLois Curfman McInnes   } else {
98228be2f97SBarry Smith     ierr = MatHeaderMerge(A,&B);CHKERRQ(ierr);
9833501a2bdSLois Curfman McInnes   }
9843a40ed3dSBarry Smith   PetscFunctionReturn(0);
985096963f5SLois Curfman McInnes }
986096963f5SLois Curfman McInnes 
98744cd7ae7SLois Curfman McInnes 
9886849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*);
989d84338a6SBarry Smith extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar);
9908965ea79SLois Curfman McInnes 
9914a2ae208SSatish Balay #undef __FUNCT__
9924994cf47SJed Brown #define __FUNCT__ "MatSetUp_MPIDense"
9934994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A)
994273d9f13SBarry Smith {
995dfbe8321SBarry Smith   PetscErrorCode ierr;
996273d9f13SBarry Smith 
997273d9f13SBarry Smith   PetscFunctionBegin;
998273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
999273d9f13SBarry Smith   PetscFunctionReturn(0);
1000273d9f13SBarry Smith }
1001273d9f13SBarry Smith 
100230716080SHong Zhang #undef __FUNCT__
1003488007eeSBarry Smith #define __FUNCT__ "MatAXPY_MPIDense"
1004488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
1005488007eeSBarry Smith {
1006488007eeSBarry Smith   PetscErrorCode ierr;
1007488007eeSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;
1008488007eeSBarry Smith 
1009488007eeSBarry Smith   PetscFunctionBegin;
1010488007eeSBarry Smith   ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr);
1011a3fa217bSJose E. Roman   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
1012488007eeSBarry Smith   PetscFunctionReturn(0);
1013488007eeSBarry Smith }
1014488007eeSBarry Smith 
1015ba337c44SJed Brown #undef __FUNCT__
1016ba337c44SJed Brown #define __FUNCT__ "MatConjugate_MPIDense"
10177087cfbeSBarry Smith PetscErrorCode  MatConjugate_MPIDense(Mat mat)
1018ba337c44SJed Brown {
1019ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)mat->data;
1020ba337c44SJed Brown   PetscErrorCode ierr;
1021ba337c44SJed Brown 
1022ba337c44SJed Brown   PetscFunctionBegin;
1023ba337c44SJed Brown   ierr = MatConjugate(a->A);CHKERRQ(ierr);
1024ba337c44SJed Brown   PetscFunctionReturn(0);
1025ba337c44SJed Brown }
1026ba337c44SJed Brown 
1027ba337c44SJed Brown #undef __FUNCT__
1028ba337c44SJed Brown #define __FUNCT__ "MatRealPart_MPIDense"
1029ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A)
1030ba337c44SJed Brown {
1031ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1032ba337c44SJed Brown   PetscErrorCode ierr;
1033ba337c44SJed Brown 
1034ba337c44SJed Brown   PetscFunctionBegin;
1035ba337c44SJed Brown   ierr = MatRealPart(a->A);CHKERRQ(ierr);
1036ba337c44SJed Brown   PetscFunctionReturn(0);
1037ba337c44SJed Brown }
1038ba337c44SJed Brown 
1039ba337c44SJed Brown #undef __FUNCT__
1040ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_MPIDense"
1041ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1042ba337c44SJed Brown {
1043ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1044ba337c44SJed Brown   PetscErrorCode ierr;
1045ba337c44SJed Brown 
1046ba337c44SJed Brown   PetscFunctionBegin;
1047ba337c44SJed Brown   ierr = MatImaginaryPart(a->A);CHKERRQ(ierr);
1048ba337c44SJed Brown   PetscFunctionReturn(0);
1049ba337c44SJed Brown }
1050ba337c44SJed Brown 
10510716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);
10520716a85fSBarry Smith #undef __FUNCT__
10530716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense"
10540716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
10550716a85fSBarry Smith {
10560716a85fSBarry Smith   PetscErrorCode ierr;
10570716a85fSBarry Smith   PetscInt       i,n;
10580716a85fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
10590716a85fSBarry Smith   PetscReal      *work;
10600716a85fSBarry Smith 
10610716a85fSBarry Smith   PetscFunctionBegin;
10620298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr);
1063785e854fSJed Brown   ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
10640716a85fSBarry Smith   ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr);
10650716a85fSBarry Smith   if (type == NORM_2) {
10660716a85fSBarry Smith     for (i=0; i<n; i++) work[i] *= work[i];
10670716a85fSBarry Smith   }
10680716a85fSBarry Smith   if (type == NORM_INFINITY) {
1069b2566f29SBarry Smith     ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr);
10700716a85fSBarry Smith   } else {
1071b2566f29SBarry Smith     ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr);
10720716a85fSBarry Smith   }
10730716a85fSBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
10740716a85fSBarry Smith   if (type == NORM_2) {
10758f1a2a5eSBarry Smith     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
10760716a85fSBarry Smith   }
10770716a85fSBarry Smith   PetscFunctionReturn(0);
10780716a85fSBarry Smith }
10790716a85fSBarry Smith 
108073a71a0fSBarry Smith #undef __FUNCT__
108173a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_MPIDense"
108273a71a0fSBarry Smith static PetscErrorCode  MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
108373a71a0fSBarry Smith {
108473a71a0fSBarry Smith   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;
108573a71a0fSBarry Smith   PetscErrorCode ierr;
108673a71a0fSBarry Smith   PetscScalar    *a;
108773a71a0fSBarry Smith   PetscInt       m,n,i;
108873a71a0fSBarry Smith 
108973a71a0fSBarry Smith   PetscFunctionBegin;
109073a71a0fSBarry Smith   ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr);
10918c778c55SBarry Smith   ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr);
109273a71a0fSBarry Smith   for (i=0; i<m*n; i++) {
109373a71a0fSBarry Smith     ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr);
109473a71a0fSBarry Smith   }
10958c778c55SBarry Smith   ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr);
109673a71a0fSBarry Smith   PetscFunctionReturn(0);
109773a71a0fSBarry Smith }
109873a71a0fSBarry Smith 
1099fd4e9aacSBarry Smith extern PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat,Mat);
1100fd4e9aacSBarry Smith 
11013b49f96aSBarry Smith #undef __FUNCT__
11023b49f96aSBarry Smith #define __FUNCT__ "MatMissingDiagonal_MPIDense"
11033b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MPIDense(Mat A,PetscBool  *missing,PetscInt *d)
11043b49f96aSBarry Smith {
11053b49f96aSBarry Smith   PetscFunctionBegin;
11063b49f96aSBarry Smith   *missing = PETSC_FALSE;
11073b49f96aSBarry Smith   PetscFunctionReturn(0);
11083b49f96aSBarry Smith }
11093b49f96aSBarry Smith 
11108965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
111109dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense,
111209dc0095SBarry Smith                                         MatGetRow_MPIDense,
111309dc0095SBarry Smith                                         MatRestoreRow_MPIDense,
111409dc0095SBarry Smith                                         MatMult_MPIDense,
111597304618SKris Buschelman                                 /*  4*/ MatMultAdd_MPIDense,
11167c922b88SBarry Smith                                         MatMultTranspose_MPIDense,
11177c922b88SBarry Smith                                         MatMultTransposeAdd_MPIDense,
11188965ea79SLois Curfman McInnes                                         0,
111909dc0095SBarry Smith                                         0,
112009dc0095SBarry Smith                                         0,
112197304618SKris Buschelman                                 /* 10*/ 0,
112209dc0095SBarry Smith                                         0,
112309dc0095SBarry Smith                                         0,
112409dc0095SBarry Smith                                         0,
112509dc0095SBarry Smith                                         MatTranspose_MPIDense,
112697304618SKris Buschelman                                 /* 15*/ MatGetInfo_MPIDense,
11276e4ee0c6SHong Zhang                                         MatEqual_MPIDense,
112809dc0095SBarry Smith                                         MatGetDiagonal_MPIDense,
11295b2fa520SLois Curfman McInnes                                         MatDiagonalScale_MPIDense,
113009dc0095SBarry Smith                                         MatNorm_MPIDense,
113197304618SKris Buschelman                                 /* 20*/ MatAssemblyBegin_MPIDense,
113209dc0095SBarry Smith                                         MatAssemblyEnd_MPIDense,
113309dc0095SBarry Smith                                         MatSetOption_MPIDense,
113409dc0095SBarry Smith                                         MatZeroEntries_MPIDense,
1135d519adbfSMatthew Knepley                                 /* 24*/ MatZeroRows_MPIDense,
1136919b68f7SBarry Smith                                         0,
113701b82886SBarry Smith                                         0,
113801b82886SBarry Smith                                         0,
113901b82886SBarry Smith                                         0,
11404994cf47SJed Brown                                 /* 29*/ MatSetUp_MPIDense,
1141273d9f13SBarry Smith                                         0,
114209dc0095SBarry Smith                                         0,
11438c778c55SBarry Smith                                         0,
11448c778c55SBarry Smith                                         0,
1145d519adbfSMatthew Knepley                                 /* 34*/ MatDuplicate_MPIDense,
114609dc0095SBarry Smith                                         0,
114709dc0095SBarry Smith                                         0,
114809dc0095SBarry Smith                                         0,
114909dc0095SBarry Smith                                         0,
1150d519adbfSMatthew Knepley                                 /* 39*/ MatAXPY_MPIDense,
11512ce60cd0SSatish Balay                                         MatGetSubMatrices_MPIDense,
115209dc0095SBarry Smith                                         0,
115309dc0095SBarry Smith                                         MatGetValues_MPIDense,
115409dc0095SBarry Smith                                         0,
1155d519adbfSMatthew Knepley                                 /* 44*/ 0,
115609dc0095SBarry Smith                                         MatScale_MPIDense,
11577d68702bSBarry Smith                                         MatShift_Basic,
115809dc0095SBarry Smith                                         0,
115909dc0095SBarry Smith                                         0,
116073a71a0fSBarry Smith                                 /* 49*/ MatSetRandom_MPIDense,
116109dc0095SBarry Smith                                         0,
116209dc0095SBarry Smith                                         0,
116309dc0095SBarry Smith                                         0,
116409dc0095SBarry Smith                                         0,
1165d519adbfSMatthew Knepley                                 /* 54*/ 0,
116609dc0095SBarry Smith                                         0,
116709dc0095SBarry Smith                                         0,
116809dc0095SBarry Smith                                         0,
116909dc0095SBarry Smith                                         0,
1170d519adbfSMatthew Knepley                                 /* 59*/ MatGetSubMatrix_MPIDense,
1171b9b97703SBarry Smith                                         MatDestroy_MPIDense,
1172b9b97703SBarry Smith                                         MatView_MPIDense,
1173357abbc8SBarry Smith                                         0,
117497304618SKris Buschelman                                         0,
1175d519adbfSMatthew Knepley                                 /* 64*/ 0,
117697304618SKris Buschelman                                         0,
117797304618SKris Buschelman                                         0,
117897304618SKris Buschelman                                         0,
117997304618SKris Buschelman                                         0,
1180d519adbfSMatthew Knepley                                 /* 69*/ 0,
118197304618SKris Buschelman                                         0,
118297304618SKris Buschelman                                         0,
118397304618SKris Buschelman                                         0,
118497304618SKris Buschelman                                         0,
1185d519adbfSMatthew Knepley                                 /* 74*/ 0,
118697304618SKris Buschelman                                         0,
118797304618SKris Buschelman                                         0,
118897304618SKris Buschelman                                         0,
118997304618SKris Buschelman                                         0,
1190d519adbfSMatthew Knepley                                 /* 79*/ 0,
119197304618SKris Buschelman                                         0,
119297304618SKris Buschelman                                         0,
119397304618SKris Buschelman                                         0,
11945bba2384SShri Abhyankar                                 /* 83*/ MatLoad_MPIDense,
1195865e5f61SKris Buschelman                                         0,
1196865e5f61SKris Buschelman                                         0,
1197865e5f61SKris Buschelman                                         0,
1198865e5f61SKris Buschelman                                         0,
1199865e5f61SKris Buschelman                                         0,
12009775878aSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
1201320f2790SHong Zhang                                 /* 89*/ MatMatMult_MPIDense_MPIDense,
1202320f2790SHong Zhang                                         MatMatMultSymbolic_MPIDense_MPIDense,
12039775878aSHong Zhang #else
12049775878aSHong Zhang                                 /* 89*/ 0,
1205865e5f61SKris Buschelman                                         0,
12069775878aSHong Zhang #endif
1207fd4e9aacSBarry Smith                                         MatMatMultNumeric_MPIDense,
12082fbe02b9SBarry Smith                                         0,
1209ba337c44SJed Brown                                         0,
1210d519adbfSMatthew Knepley                                 /* 94*/ 0,
1211865e5f61SKris Buschelman                                         0,
1212865e5f61SKris Buschelman                                         0,
1213ba337c44SJed Brown                                         0,
1214ba337c44SJed Brown                                         0,
1215ba337c44SJed Brown                                 /* 99*/ 0,
1216ba337c44SJed Brown                                         0,
1217ba337c44SJed Brown                                         0,
1218ba337c44SJed Brown                                         MatConjugate_MPIDense,
1219ba337c44SJed Brown                                         0,
1220ba337c44SJed Brown                                 /*104*/ 0,
1221ba337c44SJed Brown                                         MatRealPart_MPIDense,
1222ba337c44SJed Brown                                         MatImaginaryPart_MPIDense,
122386d161a7SShri Abhyankar                                         0,
122486d161a7SShri Abhyankar                                         0,
122586d161a7SShri Abhyankar                                 /*109*/ 0,
122686d161a7SShri Abhyankar                                         0,
122786d161a7SShri Abhyankar                                         0,
122886d161a7SShri Abhyankar                                         0,
12293b49f96aSBarry Smith                                         MatMissingDiagonal_MPIDense,
123086d161a7SShri Abhyankar                                 /*114*/ 0,
123186d161a7SShri Abhyankar                                         0,
123286d161a7SShri Abhyankar                                         0,
123386d161a7SShri Abhyankar                                         0,
123486d161a7SShri Abhyankar                                         0,
123586d161a7SShri Abhyankar                                 /*119*/ 0,
123686d161a7SShri Abhyankar                                         0,
123786d161a7SShri Abhyankar                                         0,
12380716a85fSBarry Smith                                         0,
12390716a85fSBarry Smith                                         0,
12400716a85fSBarry Smith                                 /*124*/ 0,
12413964eb88SJed Brown                                         MatGetColumnNorms_MPIDense,
12423964eb88SJed Brown                                         0,
12433964eb88SJed Brown                                         0,
12443964eb88SJed Brown                                         0,
12453964eb88SJed Brown                                 /*129*/ 0,
1246cb20be35SHong Zhang                                         MatTransposeMatMult_MPIDense_MPIDense,
1247cb20be35SHong Zhang                                         MatTransposeMatMultSymbolic_MPIDense_MPIDense,
1248cb20be35SHong Zhang                                         MatTransposeMatMultNumeric_MPIDense_MPIDense,
12493964eb88SJed Brown                                         0,
12503964eb88SJed Brown                                 /*134*/ 0,
12513964eb88SJed Brown                                         0,
12523964eb88SJed Brown                                         0,
12533964eb88SJed Brown                                         0,
12543964eb88SJed Brown                                         0,
12553964eb88SJed Brown                                 /*139*/ 0,
1256f9426fe0SMark Adams                                         0,
12573964eb88SJed Brown                                         0
1258ba337c44SJed Brown };
12598965ea79SLois Curfman McInnes 
12604a2ae208SSatish Balay #undef __FUNCT__
1261a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
12627087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1263a23d5eceSKris Buschelman {
1264a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1265dfbe8321SBarry Smith   PetscErrorCode ierr;
1266a23d5eceSKris Buschelman 
1267a23d5eceSKris Buschelman   PetscFunctionBegin;
1268a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1269a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1270a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1271a23d5eceSKris Buschelman 
1272a23d5eceSKris Buschelman   a       = (Mat_MPIDense*)mat->data;
127334ef9618SShri Abhyankar   ierr    = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
127434ef9618SShri Abhyankar   ierr    = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
127534ef9618SShri Abhyankar   a->nvec = mat->cmap->n;
127634ef9618SShri Abhyankar 
1277f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1278d0f46423SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr);
12795c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
12805c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
12813bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr);
1282a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1283a23d5eceSKris Buschelman }
1284a23d5eceSKris Buschelman 
128565b80a83SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
1286a23d5eceSKris Buschelman #undef __FUNCT__
12878baccfbdSHong Zhang #define __FUNCT__ "MatConvert_MPIDense_Elemental"
1288*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
12898baccfbdSHong Zhang {
12908ea901baSHong Zhang   Mat            mat_elemental;
12918ea901baSHong Zhang   PetscErrorCode ierr;
129232d7a744SHong Zhang   PetscScalar    *v;
129332d7a744SHong Zhang   PetscInt       m=A->rmap->n,N=A->cmap->N,rstart=A->rmap->rstart,i,*rows,*cols;
12948ea901baSHong Zhang 
12958baccfbdSHong Zhang   PetscFunctionBegin;
1296378336b6SHong Zhang   if (reuse == MAT_REUSE_MATRIX) {
1297378336b6SHong Zhang     mat_elemental = *newmat;
1298378336b6SHong Zhang     ierr = MatZeroEntries(*newmat);CHKERRQ(ierr);
1299378336b6SHong Zhang   } else {
1300378336b6SHong Zhang     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
1301378336b6SHong Zhang     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1302378336b6SHong Zhang     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
1303378336b6SHong Zhang     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
130432d7a744SHong Zhang     ierr = MatSetOption(mat_elemental,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
1305378336b6SHong Zhang   }
1306378336b6SHong Zhang 
130732d7a744SHong Zhang   ierr = PetscMalloc2(m,&rows,N,&cols);CHKERRQ(ierr);
130832d7a744SHong Zhang   for (i=0; i<N; i++) cols[i] = i;
130932d7a744SHong Zhang   for (i=0; i<m; i++) rows[i] = rstart + i;
13108ea901baSHong Zhang 
13118ea901baSHong Zhang   /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
131232d7a744SHong Zhang   ierr = MatDenseGetArray(A,&v);CHKERRQ(ierr);
131332d7a744SHong Zhang   ierr = MatSetValues(mat_elemental,m,rows,N,cols,v,ADD_VALUES);CHKERRQ(ierr);
13148ea901baSHong Zhang   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13158ea901baSHong Zhang   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
131632d7a744SHong Zhang   ierr = MatDenseRestoreArray(A,&v);CHKERRQ(ierr);
131732d7a744SHong Zhang   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
13188ea901baSHong Zhang 
1319511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
132028be2f97SBarry Smith     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
13218ea901baSHong Zhang   } else {
13228ea901baSHong Zhang     *newmat = mat_elemental;
13238ea901baSHong Zhang   }
13248baccfbdSHong Zhang   PetscFunctionReturn(0);
13258baccfbdSHong Zhang }
132665b80a83SHong Zhang #endif
13278baccfbdSHong Zhang 
13288baccfbdSHong Zhang #undef __FUNCT__
13294a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
13308cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat)
1331273d9f13SBarry Smith {
1332273d9f13SBarry Smith   Mat_MPIDense   *a;
1333dfbe8321SBarry Smith   PetscErrorCode ierr;
1334273d9f13SBarry Smith 
1335273d9f13SBarry Smith   PetscFunctionBegin;
1336b00a9115SJed Brown   ierr      = PetscNewLog(mat,&a);CHKERRQ(ierr);
1337b0a32e0cSBarry Smith   mat->data = (void*)a;
1338273d9f13SBarry Smith   ierr      = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1339273d9f13SBarry Smith 
1340273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1341ce94432eSBarry Smith   ierr            = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);CHKERRQ(ierr);
1342ce94432eSBarry Smith   ierr            = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);CHKERRQ(ierr);
1343273d9f13SBarry Smith 
1344273d9f13SBarry Smith   /* build cache for off array entries formed */
1345273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
13462205254eSKarl Rupp 
1347ce94432eSBarry Smith   ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);CHKERRQ(ierr);
1348273d9f13SBarry Smith 
1349273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1350273d9f13SBarry Smith   a->lvec        = 0;
1351273d9f13SBarry Smith   a->Mvctx       = 0;
1352273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1353273d9f13SBarry Smith 
1354bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);CHKERRQ(ierr);
1355bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr);
13568baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
13578baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",MatConvert_MPIDense_Elemental);CHKERRQ(ierr);
13588baccfbdSHong Zhang #endif
1359bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1360bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1361bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
1362bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
1363bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
13648949adfdSHong Zhang 
13658949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",MatTransposeMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
13668949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",MatTransposeMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
13678949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",MatTransposeMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
136838aed534SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr);
1369273d9f13SBarry Smith   PetscFunctionReturn(0);
1370273d9f13SBarry Smith }
1371273d9f13SBarry Smith 
1372209238afSKris Buschelman /*MC
1373002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1374209238afSKris Buschelman 
1375209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1376209238afSKris Buschelman    and MATMPIDENSE otherwise.
1377209238afSKris Buschelman 
1378209238afSKris Buschelman    Options Database Keys:
1379209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1380209238afSKris Buschelman 
1381209238afSKris Buschelman   Level: beginner
1382209238afSKris Buschelman 
138301b82886SBarry Smith 
1384209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1385209238afSKris Buschelman M*/
1386209238afSKris Buschelman 
13874a2ae208SSatish Balay #undef __FUNCT__
13884a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1389273d9f13SBarry Smith /*@C
1390273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1391273d9f13SBarry Smith 
1392273d9f13SBarry Smith    Not collective
1393273d9f13SBarry Smith 
1394273d9f13SBarry Smith    Input Parameters:
13951c4f3114SJed Brown .  B - the matrix
13960298fd71SBarry Smith -  data - optional location of matrix data.  Set data=NULL for PETSc
1397273d9f13SBarry Smith    to control all matrix memory allocation.
1398273d9f13SBarry Smith 
1399273d9f13SBarry Smith    Notes:
1400273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1401273d9f13SBarry Smith    storage by columns.
1402273d9f13SBarry Smith 
1403273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1404273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
14050298fd71SBarry Smith    set data=NULL.
1406273d9f13SBarry Smith 
1407273d9f13SBarry Smith    Level: intermediate
1408273d9f13SBarry Smith 
1409273d9f13SBarry Smith .keywords: matrix,dense, parallel
1410273d9f13SBarry Smith 
1411273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1412273d9f13SBarry Smith @*/
14131c4f3114SJed Brown PetscErrorCode  MatMPIDenseSetPreallocation(Mat B,PetscScalar *data)
1414273d9f13SBarry Smith {
14154ac538c5SBarry Smith   PetscErrorCode ierr;
1416273d9f13SBarry Smith 
1417273d9f13SBarry Smith   PetscFunctionBegin;
14181c4f3114SJed Brown   ierr = PetscTryMethod(B,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(B,data));CHKERRQ(ierr);
1419273d9f13SBarry Smith   PetscFunctionReturn(0);
1420273d9f13SBarry Smith }
1421273d9f13SBarry Smith 
14224a2ae208SSatish Balay #undef __FUNCT__
142369b1f4b7SBarry Smith #define __FUNCT__ "MatCreateDense"
14248965ea79SLois Curfman McInnes /*@C
142569b1f4b7SBarry Smith    MatCreateDense - Creates a parallel matrix in dense format.
14268965ea79SLois Curfman McInnes 
1427db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1428db81eaa0SLois Curfman McInnes 
14298965ea79SLois Curfman McInnes    Input Parameters:
1430db81eaa0SLois Curfman McInnes +  comm - MPI communicator
14318965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1432db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
14338965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1434db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
14356cfe35ddSJose E. Roman -  data - optional location of matrix data.  Set data=NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1436dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
14378965ea79SLois Curfman McInnes 
14388965ea79SLois Curfman McInnes    Output Parameter:
1439477f1c0bSLois Curfman McInnes .  A - the matrix
14408965ea79SLois Curfman McInnes 
1441b259b22eSLois Curfman McInnes    Notes:
144239ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
144339ddd567SLois Curfman McInnes    storage by columns.
14448965ea79SLois Curfman McInnes 
144518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
144618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
14476cfe35ddSJose E. Roman    set data=NULL (PETSC_NULL_SCALAR for Fortran users).
144818f449edSLois Curfman McInnes 
14498965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
14508965ea79SLois Curfman McInnes    (possibly both).
14518965ea79SLois Curfman McInnes 
1452027ccd11SLois Curfman McInnes    Level: intermediate
1453027ccd11SLois Curfman McInnes 
145439ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
14558965ea79SLois Curfman McInnes 
145639ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
14578965ea79SLois Curfman McInnes @*/
145869b1f4b7SBarry Smith PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
14598965ea79SLois Curfman McInnes {
14606849ba73SBarry Smith   PetscErrorCode ierr;
146113f74950SBarry Smith   PetscMPIInt    size;
14628965ea79SLois Curfman McInnes 
14633a40ed3dSBarry Smith   PetscFunctionBegin;
1464f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1465f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1466273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1467273d9f13SBarry Smith   if (size > 1) {
1468273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1469273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
14706cfe35ddSJose E. Roman     if (data) {  /* user provided data array, so no need to assemble */
14716cfe35ddSJose E. Roman       ierr = MatSetUpMultiply_MPIDense(*A);CHKERRQ(ierr);
14726cfe35ddSJose E. Roman       (*A)->assembled = PETSC_TRUE;
14736cfe35ddSJose E. Roman     }
1474273d9f13SBarry Smith   } else {
1475273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1476273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
14778c469469SLois Curfman McInnes   }
14783a40ed3dSBarry Smith   PetscFunctionReturn(0);
14798965ea79SLois Curfman McInnes }
14808965ea79SLois Curfman McInnes 
14814a2ae208SSatish Balay #undef __FUNCT__
14824a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
14836849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
14848965ea79SLois Curfman McInnes {
14858965ea79SLois Curfman McInnes   Mat            mat;
14863501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1487dfbe8321SBarry Smith   PetscErrorCode ierr;
14888965ea79SLois Curfman McInnes 
14893a40ed3dSBarry Smith   PetscFunctionBegin;
14908965ea79SLois Curfman McInnes   *newmat = 0;
1491ce94432eSBarry Smith   ierr    = MatCreate(PetscObjectComm((PetscObject)A),&mat);CHKERRQ(ierr);
1492d0f46423SBarry Smith   ierr    = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
14937adad957SLisandro Dalcin   ierr    = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1494834f8fabSBarry Smith   a       = (Mat_MPIDense*)mat->data;
1495e04c1aa4SHong Zhang   ierr    = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
14965aa7edbeSHong Zhang 
1497d5f3da31SBarry Smith   mat->factortype   = A->factortype;
1498c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1499273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
15008965ea79SLois Curfman McInnes 
15018965ea79SLois Curfman McInnes   a->size         = oldmat->size;
15028965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1503e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1504b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
15053782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1506e04c1aa4SHong Zhang 
15071e1e43feSBarry Smith   ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr);
15081e1e43feSBarry Smith   ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr);
15098965ea79SLois Curfman McInnes 
1510329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
15115609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
15123bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr);
151301b82886SBarry Smith 
15148965ea79SLois Curfman McInnes   *newmat = mat;
15153a40ed3dSBarry Smith   PetscFunctionReturn(0);
15168965ea79SLois Curfman McInnes }
15178965ea79SLois Curfman McInnes 
15184a2ae208SSatish Balay #undef __FUNCT__
15195bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
15209d36ed5fSBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat)
152186d161a7SShri Abhyankar {
152286d161a7SShri Abhyankar   PetscErrorCode ierr;
152386d161a7SShri Abhyankar   PetscMPIInt    rank,size;
152474c13d95SBarry Smith   const PetscInt *rowners;
152574c13d95SBarry Smith   PetscInt       i,m,n,nz,j,mMax;
152686d161a7SShri Abhyankar   PetscScalar    *array,*vals,*vals_ptr;
15279d36ed5fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)newmat->data;
152886d161a7SShri Abhyankar 
152986d161a7SShri Abhyankar   PetscFunctionBegin;
153086d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
153186d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
153286d161a7SShri Abhyankar 
153374c13d95SBarry Smith   /* determine ownership of rows and columns */
153474c13d95SBarry Smith   m = (newmat->rmap->n < 0) ? PETSC_DECIDE : newmat->rmap->n;
153574c13d95SBarry Smith   n = (newmat->cmap->n < 0) ? PETSC_DECIDE : newmat->cmap->n;
153686d161a7SShri Abhyankar 
153774c13d95SBarry Smith   ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr);
15389d36ed5fSBarry Smith   if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) {
15390298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
15409d36ed5fSBarry Smith   }
15418c778c55SBarry Smith   ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr);
154249467e7dSBarry Smith   ierr = MatGetLocalSize(newmat,&m,NULL);CHKERRQ(ierr);
154374c13d95SBarry Smith   ierr = MatGetOwnershipRanges(newmat,&rowners);CHKERRQ(ierr);
1544396e826eSBarry Smith   ierr = MPI_Reduce(&m,&mMax,1,MPIU_INT,MPI_MAX,0,comm);CHKERRQ(ierr);
154586d161a7SShri Abhyankar   if (!rank) {
15469d36ed5fSBarry Smith     ierr = PetscMalloc1(mMax*N,&vals);CHKERRQ(ierr);
154786d161a7SShri Abhyankar 
154886d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
154986d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
155086d161a7SShri Abhyankar 
155186d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
155286d161a7SShri Abhyankar     vals_ptr = vals;
155386d161a7SShri Abhyankar     for (i=0; i<m; i++) {
155486d161a7SShri Abhyankar       for (j=0; j<N; j++) {
155586d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
155686d161a7SShri Abhyankar       }
155786d161a7SShri Abhyankar     }
155886d161a7SShri Abhyankar 
155986d161a7SShri Abhyankar     /* read in other processors and ship out */
156086d161a7SShri Abhyankar     for (i=1; i<size; i++) {
156186d161a7SShri Abhyankar       nz   = (rowners[i+1] - rowners[i])*N;
156286d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1563a25532f0SBarry Smith       ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
156486d161a7SShri Abhyankar     }
156586d161a7SShri Abhyankar   } else {
156686d161a7SShri Abhyankar     /* receive numeric values */
1567785e854fSJed Brown     ierr = PetscMalloc1(m*N,&vals);CHKERRQ(ierr);
156886d161a7SShri Abhyankar 
156986d161a7SShri Abhyankar     /* receive message of values*/
1570a25532f0SBarry Smith     ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
157186d161a7SShri Abhyankar 
157286d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
157386d161a7SShri Abhyankar     vals_ptr = vals;
157486d161a7SShri Abhyankar     for (i=0; i<m; i++) {
157586d161a7SShri Abhyankar       for (j=0; j<N; j++) {
157686d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
157786d161a7SShri Abhyankar       }
157886d161a7SShri Abhyankar     }
157986d161a7SShri Abhyankar   }
15808c778c55SBarry Smith   ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr);
158186d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
158286d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158386d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158486d161a7SShri Abhyankar   PetscFunctionReturn(0);
158586d161a7SShri Abhyankar }
158686d161a7SShri Abhyankar 
158786d161a7SShri Abhyankar #undef __FUNCT__
15885bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense"
1589112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer)
159086d161a7SShri Abhyankar {
1591dfdea288SBarry Smith   Mat_MPIDense   *a;
159286d161a7SShri Abhyankar   PetscScalar    *vals,*svals;
1593ce94432eSBarry Smith   MPI_Comm       comm;
159486d161a7SShri Abhyankar   MPI_Status     status;
159549467e7dSBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,n,maxnz;
159686d161a7SShri Abhyankar   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
159749467e7dSBarry Smith   PetscInt       *ourlens,*procsnz = 0,jj,*mycols,*smycols;
15989d36ed5fSBarry Smith   PetscInt       i,nz,j,rstart,rend;
159986d161a7SShri Abhyankar   int            fd;
160086d161a7SShri Abhyankar   PetscErrorCode ierr;
160186d161a7SShri Abhyankar 
160286d161a7SShri Abhyankar   PetscFunctionBegin;
1603c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
1604c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
1605ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
160686d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
160786d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
160886d161a7SShri Abhyankar   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
16095872f025SBarry Smith   if (!rank) {
161086d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,(char*)header,4,PETSC_INT);CHKERRQ(ierr);
161186d161a7SShri Abhyankar     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
161286d161a7SShri Abhyankar   }
161386d161a7SShri Abhyankar   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
161486d161a7SShri Abhyankar   M    = header[1]; N = header[2]; nz = header[3];
161586d161a7SShri Abhyankar 
161686d161a7SShri Abhyankar   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
16179d36ed5fSBarry Smith   if (newmat->rmap->N < 0) newmat->rmap->N = M;
16189d36ed5fSBarry Smith   if (newmat->cmap->N < 0) newmat->cmap->N = N;
161986d161a7SShri Abhyankar 
16209d36ed5fSBarry Smith   if (newmat->rmap->N != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%D) and input matrix has (%D)",M,newmat->rmap->N);
16219d36ed5fSBarry Smith   if (newmat->cmap->N != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%D) and input matrix has (%D)",N,newmat->cmap->N);
162286d161a7SShri Abhyankar 
162386d161a7SShri Abhyankar   /*
162486d161a7SShri Abhyankar        Handle case where matrix is stored on disk as a dense matrix
162586d161a7SShri Abhyankar   */
162686d161a7SShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
16279d36ed5fSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
162886d161a7SShri Abhyankar     PetscFunctionReturn(0);
162986d161a7SShri Abhyankar   }
163086d161a7SShri Abhyankar 
163186d161a7SShri Abhyankar   /* determine ownership of all rows */
16322205254eSKarl Rupp   if (newmat->rmap->n < 0) {
16332205254eSKarl Rupp     ierr = PetscMPIIntCast(M/size + ((M % size) > rank),&m);CHKERRQ(ierr);
16342205254eSKarl Rupp   } else {
16352205254eSKarl Rupp     ierr = PetscMPIIntCast(newmat->rmap->n,&m);CHKERRQ(ierr);
16362205254eSKarl Rupp   }
163749467e7dSBarry Smith   if (newmat->cmap->n < 0) {
163849467e7dSBarry Smith     n = PETSC_DECIDE;
163949467e7dSBarry Smith   } else {
164049467e7dSBarry Smith     ierr = PetscMPIIntCast(newmat->cmap->n,&n);CHKERRQ(ierr);
164149467e7dSBarry Smith   }
164249467e7dSBarry Smith 
1643854ce69bSBarry Smith   ierr       = PetscMalloc1(size+2,&rowners);CHKERRQ(ierr);
164486d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
164586d161a7SShri Abhyankar   rowners[0] = 0;
164686d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
164786d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
164886d161a7SShri Abhyankar   }
164986d161a7SShri Abhyankar   rstart = rowners[rank];
165086d161a7SShri Abhyankar   rend   = rowners[rank+1];
165186d161a7SShri Abhyankar 
165286d161a7SShri Abhyankar   /* distribute row lengths to all processors */
165349467e7dSBarry Smith   ierr = PetscMalloc1(rend-rstart,&ourlens);CHKERRQ(ierr);
165486d161a7SShri Abhyankar   if (!rank) {
1655785e854fSJed Brown     ierr = PetscMalloc1(M,&rowlengths);CHKERRQ(ierr);
165686d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1657785e854fSJed Brown     ierr = PetscMalloc1(size,&sndcounts);CHKERRQ(ierr);
165886d161a7SShri Abhyankar     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
165986d161a7SShri Abhyankar     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
166086d161a7SShri Abhyankar     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
166186d161a7SShri Abhyankar   } else {
166286d161a7SShri Abhyankar     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
166386d161a7SShri Abhyankar   }
166486d161a7SShri Abhyankar 
166586d161a7SShri Abhyankar   if (!rank) {
166686d161a7SShri Abhyankar     /* calculate the number of nonzeros on each processor */
1667785e854fSJed Brown     ierr = PetscMalloc1(size,&procsnz);CHKERRQ(ierr);
166886d161a7SShri Abhyankar     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
166986d161a7SShri Abhyankar     for (i=0; i<size; i++) {
167086d161a7SShri Abhyankar       for (j=rowners[i]; j< rowners[i+1]; j++) {
167186d161a7SShri Abhyankar         procsnz[i] += rowlengths[j];
167286d161a7SShri Abhyankar       }
167386d161a7SShri Abhyankar     }
167486d161a7SShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
167586d161a7SShri Abhyankar 
167686d161a7SShri Abhyankar     /* determine max buffer needed and allocate it */
167786d161a7SShri Abhyankar     maxnz = 0;
167886d161a7SShri Abhyankar     for (i=0; i<size; i++) {
167986d161a7SShri Abhyankar       maxnz = PetscMax(maxnz,procsnz[i]);
168086d161a7SShri Abhyankar     }
1681785e854fSJed Brown     ierr = PetscMalloc1(maxnz,&cols);CHKERRQ(ierr);
168286d161a7SShri Abhyankar 
168386d161a7SShri Abhyankar     /* read in my part of the matrix column indices  */
168486d161a7SShri Abhyankar     nz   = procsnz[0];
1685785e854fSJed Brown     ierr = PetscMalloc1(nz,&mycols);CHKERRQ(ierr);
168686d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
168786d161a7SShri Abhyankar 
168886d161a7SShri Abhyankar     /* read in every one elses and ship off */
168986d161a7SShri Abhyankar     for (i=1; i<size; i++) {
169086d161a7SShri Abhyankar       nz   = procsnz[i];
169186d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
169286d161a7SShri Abhyankar       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
169386d161a7SShri Abhyankar     }
169486d161a7SShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
169586d161a7SShri Abhyankar   } else {
169686d161a7SShri Abhyankar     /* determine buffer space needed for message */
169786d161a7SShri Abhyankar     nz = 0;
169886d161a7SShri Abhyankar     for (i=0; i<m; i++) {
169986d161a7SShri Abhyankar       nz += ourlens[i];
170086d161a7SShri Abhyankar     }
1701854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&mycols);CHKERRQ(ierr);
170286d161a7SShri Abhyankar 
170386d161a7SShri Abhyankar     /* receive message of column indices*/
170486d161a7SShri Abhyankar     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
170586d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
170686d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
170786d161a7SShri Abhyankar   }
170886d161a7SShri Abhyankar 
170949467e7dSBarry Smith   ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr);
1710dfdea288SBarry Smith   a = (Mat_MPIDense*)newmat->data;
17112e3566b0SBarry Smith   if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) {
17120298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
1713dfdea288SBarry Smith   }
171486d161a7SShri Abhyankar 
171586d161a7SShri Abhyankar   if (!rank) {
1716785e854fSJed Brown     ierr = PetscMalloc1(maxnz,&vals);CHKERRQ(ierr);
171786d161a7SShri Abhyankar 
171886d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
171986d161a7SShri Abhyankar     nz   = procsnz[0];
172086d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
172186d161a7SShri Abhyankar 
172286d161a7SShri Abhyankar     /* insert into matrix */
172386d161a7SShri Abhyankar     jj      = rstart;
172486d161a7SShri Abhyankar     smycols = mycols;
172586d161a7SShri Abhyankar     svals   = vals;
172686d161a7SShri Abhyankar     for (i=0; i<m; i++) {
172786d161a7SShri Abhyankar       ierr     = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
172886d161a7SShri Abhyankar       smycols += ourlens[i];
172986d161a7SShri Abhyankar       svals   += ourlens[i];
173086d161a7SShri Abhyankar       jj++;
173186d161a7SShri Abhyankar     }
173286d161a7SShri Abhyankar 
173386d161a7SShri Abhyankar     /* read in other processors and ship out */
173486d161a7SShri Abhyankar     for (i=1; i<size; i++) {
173586d161a7SShri Abhyankar       nz   = procsnz[i];
173686d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
173786d161a7SShri Abhyankar       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr);
173886d161a7SShri Abhyankar     }
173986d161a7SShri Abhyankar     ierr = PetscFree(procsnz);CHKERRQ(ierr);
174086d161a7SShri Abhyankar   } else {
174186d161a7SShri Abhyankar     /* receive numeric values */
1742854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&vals);CHKERRQ(ierr);
174386d161a7SShri Abhyankar 
174486d161a7SShri Abhyankar     /* receive message of values*/
174586d161a7SShri Abhyankar     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr);
174686d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
174786d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
174886d161a7SShri Abhyankar 
174986d161a7SShri Abhyankar     /* insert into matrix */
175086d161a7SShri Abhyankar     jj      = rstart;
175186d161a7SShri Abhyankar     smycols = mycols;
175286d161a7SShri Abhyankar     svals   = vals;
175386d161a7SShri Abhyankar     for (i=0; i<m; i++) {
175486d161a7SShri Abhyankar       ierr     = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
175586d161a7SShri Abhyankar       smycols += ourlens[i];
175686d161a7SShri Abhyankar       svals   += ourlens[i];
175786d161a7SShri Abhyankar       jj++;
175886d161a7SShri Abhyankar     }
175986d161a7SShri Abhyankar   }
176049467e7dSBarry Smith   ierr = PetscFree(ourlens);CHKERRQ(ierr);
176186d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
176286d161a7SShri Abhyankar   ierr = PetscFree(mycols);CHKERRQ(ierr);
176386d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
176486d161a7SShri Abhyankar 
176586d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
176686d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
176786d161a7SShri Abhyankar   PetscFunctionReturn(0);
176886d161a7SShri Abhyankar }
176986d161a7SShri Abhyankar 
177086d161a7SShri Abhyankar #undef __FUNCT__
17716e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense"
1772ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
17736e4ee0c6SHong Zhang {
17746e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
17756e4ee0c6SHong Zhang   Mat            a,b;
1776ace3abfcSBarry Smith   PetscBool      flg;
17776e4ee0c6SHong Zhang   PetscErrorCode ierr;
177890ace30eSBarry Smith 
17796e4ee0c6SHong Zhang   PetscFunctionBegin;
17806e4ee0c6SHong Zhang   a    = matA->A;
17816e4ee0c6SHong Zhang   b    = matB->A;
17826e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
1783b2566f29SBarry Smith   ierr = MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
17846e4ee0c6SHong Zhang   PetscFunctionReturn(0);
17856e4ee0c6SHong Zhang }
178690ace30eSBarry Smith 
1787cb20be35SHong Zhang #undef __FUNCT__
1788baa3c1c6SHong Zhang #define __FUNCT__ "MatDestroy_MatTransMatMult_MPIDense_MPIDense"
1789baa3c1c6SHong Zhang PetscErrorCode MatDestroy_MatTransMatMult_MPIDense_MPIDense(Mat A)
1790baa3c1c6SHong Zhang {
1791baa3c1c6SHong Zhang   PetscErrorCode        ierr;
1792baa3c1c6SHong Zhang   Mat_MPIDense          *a = (Mat_MPIDense*)A->data;
1793baa3c1c6SHong Zhang   Mat_TransMatMultDense *atb = a->atbdense;
1794baa3c1c6SHong Zhang 
1795baa3c1c6SHong Zhang   PetscFunctionBegin;
1796c5ef1628SHong Zhang   ierr = PetscFree3(atb->sendbuf,atb->atbarray,atb->recvcounts);CHKERRQ(ierr);
1797baa3c1c6SHong Zhang   ierr = (atb->destroy)(A);CHKERRQ(ierr);
1798baa3c1c6SHong Zhang   ierr = PetscFree(atb);CHKERRQ(ierr);
1799baa3c1c6SHong Zhang   PetscFunctionReturn(0);
1800baa3c1c6SHong Zhang }
1801baa3c1c6SHong Zhang 
1802baa3c1c6SHong Zhang #undef __FUNCT__
1803cb20be35SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIDense_MPIDense"
1804cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C)
1805cb20be35SHong Zhang {
1806baa3c1c6SHong Zhang   Mat_MPIDense   *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data;
1807660d5466SHong Zhang   Mat_SeqDense   *aseq=(Mat_SeqDense*)(a->A)->data, *bseq=(Mat_SeqDense*)(b->A)->data;
1808baa3c1c6SHong Zhang   Mat_TransMatMultDense *atb = c->atbdense;
1809cb20be35SHong Zhang   PetscErrorCode ierr;
1810cb20be35SHong Zhang   MPI_Comm       comm;
1811d5017740SHong Zhang   PetscMPIInt    rank,size,*recvcounts=atb->recvcounts;
1812c5ef1628SHong Zhang   PetscScalar    *carray,*atbarray=atb->atbarray,*sendbuf=atb->sendbuf;
1813d5017740SHong Zhang   PetscInt       i,cN=C->cmap->N,cM=C->rmap->N,proc,k,j;
1814660d5466SHong Zhang   PetscScalar    _DOne=1.0,_DZero=0.0;
1815adc7a786SHong Zhang   PetscBLASInt   am,an,bn,aN;
1816e68c0b26SHong Zhang   const PetscInt *ranges;
1817cb20be35SHong Zhang 
1818cb20be35SHong Zhang   PetscFunctionBegin;
1819cb20be35SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1820cb20be35SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1821cb20be35SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1822e68c0b26SHong Zhang 
1823c5ef1628SHong Zhang   /* compute atbarray = aseq^T * bseq */
1824660d5466SHong Zhang   ierr = PetscBLASIntCast(a->A->cmap->n,&an);CHKERRQ(ierr);
1825660d5466SHong Zhang   ierr = PetscBLASIntCast(b->A->cmap->n,&bn);CHKERRQ(ierr);
1826660d5466SHong Zhang   ierr = PetscBLASIntCast(a->A->rmap->n,&am);CHKERRQ(ierr);
1827adc7a786SHong Zhang   ierr = PetscBLASIntCast(A->cmap->N,&aN);CHKERRQ(ierr);
1828adc7a786SHong Zhang   PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&an,&bn,&am,&_DOne,aseq->v,&aseq->lda,bseq->v,&bseq->lda,&_DZero,atbarray,&aN));
1829cb20be35SHong Zhang 
1830cb20be35SHong Zhang   ierr = MatGetOwnershipRanges(C,&ranges);CHKERRQ(ierr);
1831c5ef1628SHong Zhang   for (i=0; i<size; i++) recvcounts[i] = (ranges[i+1] - ranges[i])*cN;
1832cb20be35SHong Zhang 
1833660d5466SHong Zhang   /* arrange atbarray into sendbuf */
1834cb20be35SHong Zhang   k = 0;
1835cb20be35SHong Zhang   for (proc=0; proc<size; proc++) {
1836baa3c1c6SHong Zhang     for (j=0; j<cN; j++) {
1837c5ef1628SHong Zhang       for (i=ranges[proc]; i<ranges[proc+1]; i++) sendbuf[k++] = atbarray[i+j*cM];
1838cb20be35SHong Zhang     }
1839cb20be35SHong Zhang   }
1840c5ef1628SHong Zhang   /* sum all atbarray to local values of C */
1841660d5466SHong Zhang   ierr = MatDenseGetArray(c->A,&carray);CHKERRQ(ierr);
18423462b7efSHong Zhang   ierr = MPI_Reduce_scatter(sendbuf,carray,recvcounts,MPIU_SCALAR,MPIU_SUM,comm);CHKERRQ(ierr);
1843660d5466SHong Zhang   ierr = MatDenseRestoreArray(c->A,&carray);CHKERRQ(ierr);
1844cb20be35SHong Zhang   PetscFunctionReturn(0);
1845cb20be35SHong Zhang }
1846cb20be35SHong Zhang 
1847cb20be35SHong Zhang #undef __FUNCT__
1848cb20be35SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIDense_MPIDense"
1849cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
1850cb20be35SHong Zhang {
1851cb20be35SHong Zhang   PetscErrorCode        ierr;
1852baa3c1c6SHong Zhang   Mat                   Cdense;
1853cb20be35SHong Zhang   MPI_Comm              comm;
1854baa3c1c6SHong Zhang   PetscMPIInt           size;
1855660d5466SHong Zhang   PetscInt              cm=A->cmap->n,cM,cN=B->cmap->N;
1856baa3c1c6SHong Zhang   Mat_MPIDense          *c;
1857baa3c1c6SHong Zhang   Mat_TransMatMultDense *atb;
1858cb20be35SHong Zhang 
1859cb20be35SHong Zhang   PetscFunctionBegin;
1860baa3c1c6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1861cb20be35SHong Zhang   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) {
1862cb20be35SHong Zhang     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);
1863cb20be35SHong Zhang   }
1864cb20be35SHong Zhang 
1865baa3c1c6SHong Zhang   /* create matrix product Cdense */
1866baa3c1c6SHong Zhang   ierr = MatCreate(comm,&Cdense);CHKERRQ(ierr);
1867660d5466SHong Zhang   ierr = MatSetSizes(Cdense,cm,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1868baa3c1c6SHong Zhang   ierr = MatSetType(Cdense,MATMPIDENSE);CHKERRQ(ierr);
1869baa3c1c6SHong Zhang   ierr = MatMPIDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr);
1870baa3c1c6SHong Zhang   ierr = MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1871baa3c1c6SHong Zhang   ierr = MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1872baa3c1c6SHong Zhang   *C   = Cdense;
1873baa3c1c6SHong Zhang 
1874baa3c1c6SHong Zhang   /* create data structure for reuse Cdense */
1875baa3c1c6SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1876baa3c1c6SHong Zhang   ierr = PetscNew(&atb);CHKERRQ(ierr);
1877660d5466SHong Zhang   cM = Cdense->rmap->N;
1878c5ef1628SHong Zhang   ierr = PetscMalloc3(cM*cN,&atb->sendbuf,cM*cN,&atb->atbarray,size,&atb->recvcounts);CHKERRQ(ierr);
1879baa3c1c6SHong Zhang 
1880baa3c1c6SHong Zhang   c                    = (Mat_MPIDense*)Cdense->data;
1881baa3c1c6SHong Zhang   c->atbdense          = atb;
1882baa3c1c6SHong Zhang   atb->destroy         = Cdense->ops->destroy;
1883baa3c1c6SHong Zhang   Cdense->ops->destroy = MatDestroy_MatTransMatMult_MPIDense_MPIDense;
1884cb20be35SHong Zhang   PetscFunctionReturn(0);
1885cb20be35SHong Zhang }
1886cb20be35SHong Zhang 
1887cb20be35SHong Zhang #undef __FUNCT__
1888cb20be35SHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIDense_MPIDense"
1889cb20be35SHong Zhang PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1890cb20be35SHong Zhang {
1891cb20be35SHong Zhang   PetscErrorCode ierr;
1892cb20be35SHong Zhang 
1893cb20be35SHong Zhang   PetscFunctionBegin;
1894cb20be35SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
1895cb20be35SHong Zhang     ierr = MatTransposeMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr);
1896cb20be35SHong Zhang   }
1897cb20be35SHong Zhang   ierr = MatTransposeMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr);
1898cb20be35SHong Zhang   PetscFunctionReturn(0);
1899cb20be35SHong Zhang }
1900320f2790SHong Zhang 
1901320f2790SHong Zhang #undef __FUNCT__
1902320f2790SHong Zhang #define __FUNCT__ "MatDestroy_MatMatMult_MPIDense_MPIDense"
1903320f2790SHong Zhang PetscErrorCode MatDestroy_MatMatMult_MPIDense_MPIDense(Mat A)
1904320f2790SHong Zhang {
1905320f2790SHong Zhang   PetscErrorCode   ierr;
1906320f2790SHong Zhang   Mat_MPIDense     *a = (Mat_MPIDense*)A->data;
1907320f2790SHong Zhang   Mat_MatMultDense *ab = a->abdense;
1908320f2790SHong Zhang 
1909320f2790SHong Zhang   PetscFunctionBegin;
1910320f2790SHong Zhang   ierr = MatDestroy(&ab->Ce);CHKERRQ(ierr);
1911320f2790SHong Zhang   ierr = MatDestroy(&ab->Ae);CHKERRQ(ierr);
1912320f2790SHong Zhang   ierr = MatDestroy(&ab->Be);CHKERRQ(ierr);
1913320f2790SHong Zhang 
1914320f2790SHong Zhang   ierr = (ab->destroy)(A);CHKERRQ(ierr);
1915320f2790SHong Zhang   ierr = PetscFree(ab);CHKERRQ(ierr);
1916320f2790SHong Zhang   PetscFunctionReturn(0);
1917320f2790SHong Zhang }
1918320f2790SHong Zhang 
19195242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
1920320f2790SHong Zhang #undef __FUNCT__
1921320f2790SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIDense_MPIDense"
1922320f2790SHong Zhang PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C)
1923320f2790SHong Zhang {
1924320f2790SHong Zhang   PetscErrorCode   ierr;
1925320f2790SHong Zhang   Mat_MPIDense     *c=(Mat_MPIDense*)C->data;
1926320f2790SHong Zhang   Mat_MatMultDense *ab=c->abdense;
1927320f2790SHong Zhang 
1928320f2790SHong Zhang   PetscFunctionBegin;
1929de0a22f0SHong Zhang   ierr = MatConvert_MPIDense_Elemental(A,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Ae);CHKERRQ(ierr);
1930de0a22f0SHong Zhang   ierr = MatConvert_MPIDense_Elemental(B,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Be);CHKERRQ(ierr);
1931320f2790SHong Zhang   ierr = MatMatMultNumeric(ab->Ae,ab->Be,ab->Ce);CHKERRQ(ierr);
1932de0a22f0SHong Zhang   ierr = MatConvert(ab->Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
1933320f2790SHong Zhang   PetscFunctionReturn(0);
1934320f2790SHong Zhang }
1935320f2790SHong Zhang 
1936320f2790SHong Zhang #undef __FUNCT__
1937320f2790SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIDense_MPIDense"
1938320f2790SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
1939320f2790SHong Zhang {
1940320f2790SHong Zhang   PetscErrorCode   ierr;
1941320f2790SHong Zhang   Mat              Ae,Be,Ce;
1942320f2790SHong Zhang   Mat_MPIDense     *c;
1943320f2790SHong Zhang   Mat_MatMultDense *ab;
1944320f2790SHong Zhang 
1945320f2790SHong Zhang   PetscFunctionBegin;
1946320f2790SHong Zhang   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) {
1947378336b6SHong Zhang     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);
1948320f2790SHong Zhang   }
1949320f2790SHong Zhang 
1950320f2790SHong Zhang   /* convert A and B to Elemental matrices Ae and Be */
1951320f2790SHong Zhang   ierr = MatConvert(A,MATELEMENTAL,MAT_INITIAL_MATRIX, &Ae);CHKERRQ(ierr);
1952320f2790SHong Zhang   ierr = MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX, &Be);CHKERRQ(ierr);
1953320f2790SHong Zhang 
1954320f2790SHong Zhang   /* Ce = Ae*Be */
1955320f2790SHong Zhang   ierr = MatMatMultSymbolic(Ae,Be,fill,&Ce);CHKERRQ(ierr);
1956320f2790SHong Zhang   ierr = MatMatMultNumeric(Ae,Be,Ce);CHKERRQ(ierr);
1957320f2790SHong Zhang 
1958320f2790SHong Zhang   /* convert Ce to C */
1959320f2790SHong Zhang   ierr = MatConvert(Ce,MATMPIDENSE,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
1960320f2790SHong Zhang 
1961320f2790SHong Zhang   /* create data structure for reuse Cdense */
1962320f2790SHong Zhang   ierr = PetscNew(&ab);CHKERRQ(ierr);
1963320f2790SHong Zhang   c                  = (Mat_MPIDense*)(*C)->data;
1964320f2790SHong Zhang   c->abdense         = ab;
1965320f2790SHong Zhang 
1966320f2790SHong Zhang   ab->Ae             = Ae;
1967320f2790SHong Zhang   ab->Be             = Be;
1968320f2790SHong Zhang   ab->Ce             = Ce;
1969320f2790SHong Zhang   ab->destroy        = (*C)->ops->destroy;
1970320f2790SHong Zhang   (*C)->ops->destroy        = MatDestroy_MatMatMult_MPIDense_MPIDense;
19719775878aSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIDense_MPIDense;
1972320f2790SHong Zhang   PetscFunctionReturn(0);
1973320f2790SHong Zhang }
1974320f2790SHong Zhang 
1975320f2790SHong Zhang #undef __FUNCT__
1976320f2790SHong Zhang #define __FUNCT__ "MatMatMult_MPIDense_MPIDense"
1977320f2790SHong Zhang PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1978320f2790SHong Zhang {
1979320f2790SHong Zhang   PetscErrorCode ierr;
1980320f2790SHong Zhang 
1981320f2790SHong Zhang   PetscFunctionBegin;
1982320f2790SHong Zhang   if (scall == MAT_INITIAL_MATRIX) { /* simbolic product includes numeric product */
198357299f2fSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
1984320f2790SHong Zhang     ierr = MatMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr);
198557299f2fSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
1986320f2790SHong Zhang   } else {
198757299f2fSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
1988320f2790SHong Zhang     ierr = MatMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr);
198957299f2fSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
1990320f2790SHong Zhang   }
1991320f2790SHong Zhang   PetscFunctionReturn(0);
1992320f2790SHong Zhang }
19935242a7b1SHong Zhang #endif
1994