xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision b2566f29c2b6470df769aa9f7deb9e2726b0959e)
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>
98965ea79SLois Curfman McInnes 
10ba8c8a56SBarry Smith #undef __FUNCT__
11ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
12ab92ecdeSBarry Smith /*@
13ab92ecdeSBarry Smith 
14ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
15ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
16ab92ecdeSBarry Smith 
17ab92ecdeSBarry Smith     Input Parameter:
18ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
19ab92ecdeSBarry Smith 
20ab92ecdeSBarry Smith     Output Parameter:
21ab92ecdeSBarry Smith .      B - the inner matrix
22ab92ecdeSBarry Smith 
238e6c10adSSatish Balay     Level: intermediate
248e6c10adSSatish Balay 
25ab92ecdeSBarry Smith @*/
26ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
27ab92ecdeSBarry Smith {
28ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
29ab92ecdeSBarry Smith   PetscErrorCode ierr;
30ace3abfcSBarry Smith   PetscBool      flg;
31ab92ecdeSBarry Smith 
32ab92ecdeSBarry Smith   PetscFunctionBegin;
33251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
342205254eSKarl Rupp   if (flg) *B = mat->A;
352205254eSKarl Rupp   else *B = A;
36ab92ecdeSBarry Smith   PetscFunctionReturn(0);
37ab92ecdeSBarry Smith }
38ab92ecdeSBarry Smith 
39ab92ecdeSBarry Smith #undef __FUNCT__
40ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
41ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
42ba8c8a56SBarry Smith {
43ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
44ba8c8a56SBarry Smith   PetscErrorCode ierr;
45d0f46423SBarry Smith   PetscInt       lrow,rstart = A->rmap->rstart,rend = A->rmap->rend;
46ba8c8a56SBarry Smith 
47ba8c8a56SBarry Smith   PetscFunctionBegin;
48e7e72b3dSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows");
49ba8c8a56SBarry Smith   lrow = row - rstart;
50ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt**)idx,(const PetscScalar**)v);CHKERRQ(ierr);
51ba8c8a56SBarry Smith   PetscFunctionReturn(0);
52ba8c8a56SBarry Smith }
53ba8c8a56SBarry Smith 
54ba8c8a56SBarry Smith #undef __FUNCT__
55ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
56ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
57ba8c8a56SBarry Smith {
58ba8c8a56SBarry Smith   PetscErrorCode ierr;
59ba8c8a56SBarry Smith 
60ba8c8a56SBarry Smith   PetscFunctionBegin;
61ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
62ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
63ba8c8a56SBarry Smith   PetscFunctionReturn(0);
64ba8c8a56SBarry Smith }
65ba8c8a56SBarry Smith 
664a2ae208SSatish Balay #undef __FUNCT__
674a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
6811bd1e4dSLisandro Dalcin PetscErrorCode  MatGetDiagonalBlock_MPIDense(Mat A,Mat *a)
690de54da6SSatish Balay {
700de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
716849ba73SBarry Smith   PetscErrorCode ierr;
72d0f46423SBarry Smith   PetscInt       m = A->rmap->n,rstart = A->rmap->rstart;
7387828ca2SBarry Smith   PetscScalar    *array;
740de54da6SSatish Balay   MPI_Comm       comm;
7511bd1e4dSLisandro Dalcin   Mat            B;
760de54da6SSatish Balay 
770de54da6SSatish Balay   PetscFunctionBegin;
78e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only square matrices supported.");
790de54da6SSatish Balay 
8011bd1e4dSLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);CHKERRQ(ierr);
8111bd1e4dSLisandro Dalcin   if (!B) {
820de54da6SSatish Balay     ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
8311bd1e4dSLisandro Dalcin     ierr = MatCreate(comm,&B);CHKERRQ(ierr);
8411bd1e4dSLisandro Dalcin     ierr = MatSetSizes(B,m,m,m,m);CHKERRQ(ierr);
8511bd1e4dSLisandro Dalcin     ierr = MatSetType(B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr);
868c778c55SBarry Smith     ierr = MatDenseGetArray(mdn->A,&array);CHKERRQ(ierr);
8711bd1e4dSLisandro Dalcin     ierr = MatSeqDenseSetPreallocation(B,array+m*rstart);CHKERRQ(ierr);
888c778c55SBarry Smith     ierr = MatDenseRestoreArray(mdn->A,&array);CHKERRQ(ierr);
8911bd1e4dSLisandro Dalcin     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9011bd1e4dSLisandro Dalcin     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9111bd1e4dSLisandro Dalcin     ierr = PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);CHKERRQ(ierr);
9211bd1e4dSLisandro Dalcin     *a   = B;
9311bd1e4dSLisandro Dalcin     ierr = MatDestroy(&B);CHKERRQ(ierr);
942205254eSKarl Rupp   } else *a = B;
950de54da6SSatish Balay   PetscFunctionReturn(0);
960de54da6SSatish Balay }
970de54da6SSatish Balay 
984a2ae208SSatish Balay #undef __FUNCT__
994a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
10013f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1018965ea79SLois Curfman McInnes {
10239b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
103dfbe8321SBarry Smith   PetscErrorCode ierr;
104d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
105ace3abfcSBarry Smith   PetscBool      roworiented = A->roworiented;
1068965ea79SLois Curfman McInnes 
1073a40ed3dSBarry Smith   PetscFunctionBegin;
1088965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1095ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
110e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1118965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1128965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
11339b7565bSBarry Smith       if (roworiented) {
11439b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1153a40ed3dSBarry Smith       } else {
1168965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1175ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
118e32f2f54SBarry Smith           if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
11939b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
12039b7565bSBarry Smith         }
1218965ea79SLois Curfman McInnes       }
1222205254eSKarl Rupp     } else if (!A->donotstash) {
1235080c13bSMatthew G Knepley       mat->assembled = PETSC_FALSE;
12439b7565bSBarry Smith       if (roworiented) {
125b400d20cSBarry Smith         ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n,PETSC_FALSE);CHKERRQ(ierr);
126d36fbae8SSatish Balay       } else {
127b400d20cSBarry Smith         ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m,PETSC_FALSE);CHKERRQ(ierr);
12839b7565bSBarry Smith       }
129b49de8d1SLois Curfman McInnes     }
130b49de8d1SLois Curfman McInnes   }
1313a40ed3dSBarry Smith   PetscFunctionReturn(0);
132b49de8d1SLois Curfman McInnes }
133b49de8d1SLois Curfman McInnes 
1344a2ae208SSatish Balay #undef __FUNCT__
1354a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
13613f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
137b49de8d1SLois Curfman McInnes {
138b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
139dfbe8321SBarry Smith   PetscErrorCode ierr;
140d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
141b49de8d1SLois Curfman McInnes 
1423a40ed3dSBarry Smith   PetscFunctionBegin;
143b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
144e32f2f54SBarry Smith     if (idxm[i] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
145e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
146b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
147b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
148b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
149e32f2f54SBarry Smith         if (idxn[j] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
150e7e72b3dSBarry Smith         if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
151b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
152b49de8d1SLois Curfman McInnes       }
153e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
1548965ea79SLois Curfman McInnes   }
1553a40ed3dSBarry Smith   PetscFunctionReturn(0);
1568965ea79SLois Curfman McInnes }
1578965ea79SLois Curfman McInnes 
1584a2ae208SSatish Balay #undef __FUNCT__
1598c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_MPIDense"
1608c778c55SBarry Smith PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[])
161ff14e315SSatish Balay {
162ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
163dfbe8321SBarry Smith   PetscErrorCode ierr;
164ff14e315SSatish Balay 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
1668c778c55SBarry Smith   ierr = MatDenseGetArray(a->A,array);CHKERRQ(ierr);
1673a40ed3dSBarry Smith   PetscFunctionReturn(0);
168ff14e315SSatish Balay }
169ff14e315SSatish Balay 
1704a2ae208SSatish Balay #undef __FUNCT__
1714a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
1724aa3045dSJed Brown static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
173ca3fa75bSLois Curfman McInnes {
174ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat  = (Mat_MPIDense*)A->data,*newmatd;
175ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
1766849ba73SBarry Smith   PetscErrorCode ierr;
1774aa3045dSJed Brown   PetscInt       i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols;
1785d0c19d7SBarry Smith   const PetscInt *irow,*icol;
17987828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
180ca3fa75bSLois Curfman McInnes   Mat            newmat;
1814aa3045dSJed Brown   IS             iscol_local;
182ca3fa75bSLois Curfman McInnes 
183ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
1844aa3045dSJed Brown   ierr = ISAllGather(iscol,&iscol_local);CHKERRQ(ierr);
185ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
1864aa3045dSJed Brown   ierr = ISGetIndices(iscol_local,&icol);CHKERRQ(ierr);
187b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
188b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
1894aa3045dSJed Brown   ierr = ISGetSize(iscol,&Ncols);CHKERRQ(ierr); /* global number of columns, size of iscol_local */
190ca3fa75bSLois Curfman McInnes 
191ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1927eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1937eba5e9cSLois Curfman McInnes      original matrix! */
194ca3fa75bSLois Curfman McInnes 
195ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1967eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
197ca3fa75bSLois Curfman McInnes 
198ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
199ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
200e32f2f54SBarry Smith     /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2017eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
202ca3fa75bSLois Curfman McInnes     newmat = *B;
203ca3fa75bSLois Curfman McInnes   } else {
204ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
205ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr);
2064aa3045dSJed Brown     ierr = MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);CHKERRQ(ierr);
2077adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
2080298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
209ca3fa75bSLois Curfman McInnes   }
210ca3fa75bSLois Curfman McInnes 
211ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
212ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
213ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense*)newmatd->A->data)->v;
214ca3fa75bSLois Curfman McInnes 
2154aa3045dSJed Brown   for (i=0; i<Ncols; i++) {
21625a33276SHong Zhang     av = v + ((Mat_SeqDense*)mat->A->data)->lda*icol[i];
217ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2187eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
219ca3fa75bSLois Curfman McInnes     }
220ca3fa75bSLois Curfman McInnes   }
221ca3fa75bSLois Curfman McInnes 
222ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
223ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
224ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
225ca3fa75bSLois Curfman McInnes 
226ca3fa75bSLois Curfman McInnes   /* Free work space */
227ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
2285bdf786aSShri Abhyankar   ierr = ISRestoreIndices(iscol_local,&icol);CHKERRQ(ierr);
22932bb1f2dSLisandro Dalcin   ierr = ISDestroy(&iscol_local);CHKERRQ(ierr);
230ca3fa75bSLois Curfman McInnes   *B   = newmat;
231ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
232ca3fa75bSLois Curfman McInnes }
233ca3fa75bSLois Curfman McInnes 
2344a2ae208SSatish Balay #undef __FUNCT__
2358c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_MPIDense"
2368c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[])
237ff14e315SSatish Balay {
23873a71a0fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
23973a71a0fSBarry Smith   PetscErrorCode ierr;
24073a71a0fSBarry Smith 
2413a40ed3dSBarry Smith   PetscFunctionBegin;
2428c778c55SBarry Smith   ierr = MatDenseRestoreArray(a->A,array);CHKERRQ(ierr);
2433a40ed3dSBarry Smith   PetscFunctionReturn(0);
244ff14e315SSatish Balay }
245ff14e315SSatish Balay 
2464a2ae208SSatish Balay #undef __FUNCT__
2474a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
248dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2498965ea79SLois Curfman McInnes {
25039ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
251ce94432eSBarry Smith   MPI_Comm       comm;
252dfbe8321SBarry Smith   PetscErrorCode ierr;
25313f74950SBarry Smith   PetscInt       nstash,reallocs;
2548965ea79SLois Curfman McInnes   InsertMode     addv;
2558965ea79SLois Curfman McInnes 
2563a40ed3dSBarry Smith   PetscFunctionBegin;
257ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
2588965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
259*b2566f29SBarry Smith   ierr = MPIU_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);CHKERRQ(ierr);
260ce94432eSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
261e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2628965ea79SLois Curfman McInnes 
263d0f46423SBarry Smith   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr);
2648798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
265ae15b995SBarry Smith   ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
2663a40ed3dSBarry Smith   PetscFunctionReturn(0);
2678965ea79SLois Curfman McInnes }
2688965ea79SLois Curfman McInnes 
2694a2ae208SSatish Balay #undef __FUNCT__
2704a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
271dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2728965ea79SLois Curfman McInnes {
27339ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn=(Mat_MPIDense*)mat->data;
2746849ba73SBarry Smith   PetscErrorCode ierr;
27513f74950SBarry Smith   PetscInt       i,*row,*col,flg,j,rstart,ncols;
27613f74950SBarry Smith   PetscMPIInt    n;
27787828ca2SBarry Smith   PetscScalar    *val;
2788965ea79SLois Curfman McInnes 
2793a40ed3dSBarry Smith   PetscFunctionBegin;
2808965ea79SLois Curfman McInnes   /*  wait on receives */
2817ef1d9bdSSatish Balay   while (1) {
2828798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2837ef1d9bdSSatish Balay     if (!flg) break;
2848965ea79SLois Curfman McInnes 
2857ef1d9bdSSatish Balay     for (i=0; i<n;) {
2867ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2872205254eSKarl Rupp       for (j=i,rstart=row[j]; j<n; j++) {
2882205254eSKarl Rupp         if (row[j] != rstart) break;
2892205254eSKarl Rupp       }
2907ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2917ef1d9bdSSatish Balay       else       ncols = n-i;
2927ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2934b4eb8d3SJed Brown       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);CHKERRQ(ierr);
2947ef1d9bdSSatish Balay       i    = j;
2958965ea79SLois Curfman McInnes     }
2967ef1d9bdSSatish Balay   }
2978798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2988965ea79SLois Curfman McInnes 
29939ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
30039ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3018965ea79SLois Curfman McInnes 
3026d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
30339ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3048965ea79SLois Curfman McInnes   }
3053a40ed3dSBarry Smith   PetscFunctionReturn(0);
3068965ea79SLois Curfman McInnes }
3078965ea79SLois Curfman McInnes 
3084a2ae208SSatish Balay #undef __FUNCT__
3094a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
310dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3118965ea79SLois Curfman McInnes {
312dfbe8321SBarry Smith   PetscErrorCode ierr;
31339ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3143a40ed3dSBarry Smith 
3153a40ed3dSBarry Smith   PetscFunctionBegin;
3163a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3173a40ed3dSBarry Smith   PetscFunctionReturn(0);
3188965ea79SLois Curfman McInnes }
3198965ea79SLois Curfman McInnes 
3208965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3218965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3228965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3233501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3248965ea79SLois Curfman McInnes    routine.
3258965ea79SLois Curfman McInnes */
3264a2ae208SSatish Balay #undef __FUNCT__
3274a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
3282b40b63fSBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
3298965ea79SLois Curfman McInnes {
33039ddd567SLois Curfman McInnes   Mat_MPIDense      *l = (Mat_MPIDense*)A->data;
3316849ba73SBarry Smith   PetscErrorCode    ierr;
332d0f46423SBarry Smith   PetscInt          i,*owners = A->rmap->range;
33376ec1555SBarry Smith   PetscInt          *sizes,j,idx,nsends;
33413f74950SBarry Smith   PetscInt          nmax,*svalues,*starts,*owner,nrecvs;
3357adad957SLisandro Dalcin   PetscInt          *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source;
33613f74950SBarry Smith   PetscInt          *lens,*lrows,*values;
33713f74950SBarry Smith   PetscMPIInt       n,imdex,rank = l->rank,size = l->size;
338ce94432eSBarry Smith   MPI_Comm          comm;
3398965ea79SLois Curfman McInnes   MPI_Request       *send_waits,*recv_waits;
3408965ea79SLois Curfman McInnes   MPI_Status        recv_status,*send_status;
341ace3abfcSBarry Smith   PetscBool         found;
34297b48c8fSBarry Smith   const PetscScalar *xx;
34397b48c8fSBarry Smith   PetscScalar       *bb;
3448965ea79SLois Curfman McInnes 
3453a40ed3dSBarry Smith   PetscFunctionBegin;
346ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
347ce94432eSBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(comm,PETSC_ERR_SUP,"Only handles square matrices");
348b9679d65SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership");
3498965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
350f628708eSJed Brown   ierr = PetscCalloc1(2*size,&sizes);CHKERRQ(ierr);
351037dbc42SBarry Smith   ierr = PetscMalloc1(N+1,&owner);CHKERRQ(ierr);  /* see note*/
3528965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3538965ea79SLois Curfman McInnes     idx   = rows[i];
35435d8aa7fSBarry Smith     found = PETSC_FALSE;
3558965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3568965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
35776ec1555SBarry Smith         sizes[2*j]++; sizes[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3588965ea79SLois Curfman McInnes       }
3598965ea79SLois Curfman McInnes     }
360e32f2f54SBarry Smith     if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3618965ea79SLois Curfman McInnes   }
3622205254eSKarl Rupp   nsends = 0;
36376ec1555SBarry Smith   for (i=0; i<size; i++) nsends += sizes[2*i+1];
3648965ea79SLois Curfman McInnes 
3658965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
36676ec1555SBarry Smith   ierr = PetscMaxSum(comm,sizes,&nmax,&nrecvs);CHKERRQ(ierr);
3678965ea79SLois Curfman McInnes 
3688965ea79SLois Curfman McInnes   /* post receives:   */
369785e854fSJed Brown   ierr = PetscMalloc1((nrecvs+1)*(nmax+1),&rvalues);CHKERRQ(ierr);
370854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
3718965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
37213f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3738965ea79SLois Curfman McInnes   }
3748965ea79SLois Curfman McInnes 
3758965ea79SLois Curfman McInnes   /* do sends:
3768965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3778965ea79SLois Curfman McInnes          the ith processor
3788965ea79SLois Curfman McInnes   */
379854ce69bSBarry Smith   ierr = PetscMalloc1(N+1,&svalues);CHKERRQ(ierr);
380854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
381854ce69bSBarry Smith   ierr = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
3828965ea79SLois Curfman McInnes 
3838965ea79SLois Curfman McInnes   starts[0] = 0;
38476ec1555SBarry Smith   for (i=1; i<size; i++) starts[i] = starts[i-1] + sizes[2*i-2];
3852205254eSKarl Rupp   for (i=0; i<N; i++) svalues[starts[owner[i]]++] = rows[i];
3862205254eSKarl Rupp 
3872205254eSKarl Rupp   starts[0] = 0;
38876ec1555SBarry Smith   for (i=1; i<size+1; i++) starts[i] = starts[i-1] + sizes[2*i-2];
3898965ea79SLois Curfman McInnes   count = 0;
3908965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
39176ec1555SBarry Smith     if (sizes[2*i+1]) {
39276ec1555SBarry Smith       ierr = MPI_Isend(svalues+starts[i],sizes[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3938965ea79SLois Curfman McInnes     }
3948965ea79SLois Curfman McInnes   }
395606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3968965ea79SLois Curfman McInnes 
3978965ea79SLois Curfman McInnes   base = owners[rank];
3988965ea79SLois Curfman McInnes 
3998965ea79SLois Curfman McInnes   /*  wait on receives */
400dcca6d9dSJed Brown   ierr  = PetscMalloc2(nrecvs,&lens,nrecvs,&source);CHKERRQ(ierr);
40174ed9c26SBarry Smith   count = nrecvs;
40274ed9c26SBarry Smith   slen  = 0;
4038965ea79SLois Curfman McInnes   while (count) {
404ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4058965ea79SLois Curfman McInnes     /* unpack receives into our local space */
40613f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4072205254eSKarl Rupp 
4088965ea79SLois Curfman McInnes     source[imdex] = recv_status.MPI_SOURCE;
4098965ea79SLois Curfman McInnes     lens[imdex]   = n;
4108965ea79SLois Curfman McInnes     slen += n;
4118965ea79SLois Curfman McInnes     count--;
4128965ea79SLois Curfman McInnes   }
413606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4148965ea79SLois Curfman McInnes 
4158965ea79SLois Curfman McInnes   /* move the data into the send scatter */
416854ce69bSBarry Smith   ierr  = PetscMalloc1(slen+1,&lrows);CHKERRQ(ierr);
4178965ea79SLois Curfman McInnes   count = 0;
4188965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4198965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4208965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4218965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4228965ea79SLois Curfman McInnes     }
4238965ea79SLois Curfman McInnes   }
424606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
42574ed9c26SBarry Smith   ierr = PetscFree2(lens,source);CHKERRQ(ierr);
426606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
42776ec1555SBarry Smith   ierr = PetscFree(sizes);CHKERRQ(ierr);
4288965ea79SLois Curfman McInnes 
42997b48c8fSBarry Smith   /* fix right hand side if needed */
43097b48c8fSBarry Smith   if (x && b) {
43197b48c8fSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
43297b48c8fSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
43397b48c8fSBarry Smith     for (i=0; i<slen; i++) {
43497b48c8fSBarry Smith       bb[lrows[i]] = diag*xx[lrows[i]];
43597b48c8fSBarry Smith     }
43697b48c8fSBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
43797b48c8fSBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
43897b48c8fSBarry Smith   }
43997b48c8fSBarry Smith 
4408965ea79SLois Curfman McInnes   /* actually zap the local rows */
441b9679d65SBarry Smith   ierr = MatZeroRows(l->A,slen,lrows,0.0,0,0);CHKERRQ(ierr);
442e2eb51b1SBarry Smith   if (diag != 0.0) {
443b9679d65SBarry Smith     Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data;
444b9679d65SBarry Smith     PetscInt     m   = ll->lda, i;
445b9679d65SBarry Smith 
446b9679d65SBarry Smith     for (i=0; i<slen; i++) {
447b9679d65SBarry Smith       ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag;
448b9679d65SBarry Smith     }
449b9679d65SBarry Smith   }
450606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4518965ea79SLois Curfman McInnes 
4528965ea79SLois Curfman McInnes   /* wait on sends */
4538965ea79SLois Curfman McInnes   if (nsends) {
454785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
455ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
456606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4578965ea79SLois Curfman McInnes   }
458606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
459606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4603a40ed3dSBarry Smith   PetscFunctionReturn(0);
4618965ea79SLois Curfman McInnes }
4628965ea79SLois Curfman McInnes 
4634a2ae208SSatish Balay #undef __FUNCT__
4644a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
465dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4668965ea79SLois Curfman McInnes {
46739ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
468dfbe8321SBarry Smith   PetscErrorCode ierr;
469c456f294SBarry Smith 
4703a40ed3dSBarry Smith   PetscFunctionBegin;
471ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
472ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
47344cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4743a40ed3dSBarry Smith   PetscFunctionReturn(0);
4758965ea79SLois Curfman McInnes }
4768965ea79SLois Curfman McInnes 
4774a2ae208SSatish Balay #undef __FUNCT__
4784a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
479dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4808965ea79SLois Curfman McInnes {
48139ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
482dfbe8321SBarry Smith   PetscErrorCode ierr;
483c456f294SBarry Smith 
4843a40ed3dSBarry Smith   PetscFunctionBegin;
485ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
486ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
48744cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4883a40ed3dSBarry Smith   PetscFunctionReturn(0);
4898965ea79SLois Curfman McInnes }
4908965ea79SLois Curfman McInnes 
4914a2ae208SSatish Balay #undef __FUNCT__
4924a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
493dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
494096963f5SLois Curfman McInnes {
495096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
496dfbe8321SBarry Smith   PetscErrorCode ierr;
49787828ca2SBarry Smith   PetscScalar    zero = 0.0;
498096963f5SLois Curfman McInnes 
4993a40ed3dSBarry Smith   PetscFunctionBegin;
5002dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
5017c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
502ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
503ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5043a40ed3dSBarry Smith   PetscFunctionReturn(0);
505096963f5SLois Curfman McInnes }
506096963f5SLois Curfman McInnes 
5074a2ae208SSatish Balay #undef __FUNCT__
5084a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
509dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
510096963f5SLois Curfman McInnes {
511096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
512dfbe8321SBarry Smith   PetscErrorCode ierr;
513096963f5SLois Curfman McInnes 
5143a40ed3dSBarry Smith   PetscFunctionBegin;
5153501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
5167c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
517ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
518ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5193a40ed3dSBarry Smith   PetscFunctionReturn(0);
520096963f5SLois Curfman McInnes }
521096963f5SLois Curfman McInnes 
5224a2ae208SSatish Balay #undef __FUNCT__
5234a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
524dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5258965ea79SLois Curfman McInnes {
52639ddd567SLois Curfman McInnes   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
527096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
528dfbe8321SBarry Smith   PetscErrorCode ierr;
529d0f46423SBarry Smith   PetscInt       len,i,n,m = A->rmap->n,radd;
53087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
531ed3cc1f0SBarry Smith 
5323a40ed3dSBarry Smith   PetscFunctionBegin;
5332dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5341ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
535096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
536e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
537d0f46423SBarry Smith   len  = PetscMin(a->A->rmap->n,a->A->cmap->n);
538d0f46423SBarry Smith   radd = A->rmap->rstart*m;
53944cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
540096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
541096963f5SLois Curfman McInnes   }
5421ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5433a40ed3dSBarry Smith   PetscFunctionReturn(0);
5448965ea79SLois Curfman McInnes }
5458965ea79SLois Curfman McInnes 
5464a2ae208SSatish Balay #undef __FUNCT__
5474a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
548dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5498965ea79SLois Curfman McInnes {
5503501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
551dfbe8321SBarry Smith   PetscErrorCode ierr;
552ed3cc1f0SBarry Smith 
5533a40ed3dSBarry Smith   PetscFunctionBegin;
554aa482453SBarry Smith #if defined(PETSC_USE_LOG)
555d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N);
5568965ea79SLois Curfman McInnes #endif
5578798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
5586bf464f9SBarry Smith   ierr = MatDestroy(&mdn->A);CHKERRQ(ierr);
5596bf464f9SBarry Smith   ierr = VecDestroy(&mdn->lvec);CHKERRQ(ierr);
5606bf464f9SBarry Smith   ierr = VecScatterDestroy(&mdn->Mvctx);CHKERRQ(ierr);
56101b82886SBarry Smith 
562bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
563dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
5648baccfbdSHong Zhang 
5658baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);CHKERRQ(ierr);
5668baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);CHKERRQ(ierr);
5678baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
5688baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",NULL);CHKERRQ(ierr);
5698baccfbdSHong Zhang #endif
570bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",NULL);CHKERRQ(ierr);
571bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",NULL);CHKERRQ(ierr);
572bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
573bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
574bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5758baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5768baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5778baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5783a40ed3dSBarry Smith   PetscFunctionReturn(0);
5798965ea79SLois Curfman McInnes }
58039ddd567SLois Curfman McInnes 
5814a2ae208SSatish Balay #undef __FUNCT__
5824a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5836849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5848965ea79SLois Curfman McInnes {
58539ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
586dfbe8321SBarry Smith   PetscErrorCode    ierr;
587aa05aa95SBarry Smith   PetscViewerFormat format;
588aa05aa95SBarry Smith   int               fd;
589d0f46423SBarry Smith   PetscInt          header[4],mmax,N = mat->cmap->N,i,j,m,k;
590aa05aa95SBarry Smith   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
591578230a0SSatish Balay   PetscScalar       *work,*v,*vv;
592aa05aa95SBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;
5937056b6fcSBarry Smith 
5943a40ed3dSBarry Smith   PetscFunctionBegin;
59539ddd567SLois Curfman McInnes   if (mdn->size == 1) {
59639ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
597aa05aa95SBarry Smith   } else {
598aa05aa95SBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
599ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr);
600ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
601aa05aa95SBarry Smith 
602aa05aa95SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
603f4403165SShri Abhyankar     if (format == PETSC_VIEWER_NATIVE) {
604aa05aa95SBarry Smith 
605aa05aa95SBarry Smith       if (!rank) {
606aa05aa95SBarry Smith         /* store the matrix as a dense matrix */
6070700a824SBarry Smith         header[0] = MAT_FILE_CLASSID;
608d0f46423SBarry Smith         header[1] = mat->rmap->N;
609aa05aa95SBarry Smith         header[2] = N;
610aa05aa95SBarry Smith         header[3] = MATRIX_BINARY_FORMAT_DENSE;
611aa05aa95SBarry Smith         ierr      = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
612aa05aa95SBarry Smith 
613aa05aa95SBarry Smith         /* get largest work array needed for transposing array */
614d0f46423SBarry Smith         mmax = mat->rmap->n;
615aa05aa95SBarry Smith         for (i=1; i<size; i++) {
616d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
6178965ea79SLois Curfman McInnes         }
618785e854fSJed Brown         ierr = PetscMalloc1(mmax*N,&work);CHKERRQ(ierr);
619aa05aa95SBarry Smith 
620aa05aa95SBarry Smith         /* write out local array, by rows */
621d0f46423SBarry Smith         m = mat->rmap->n;
622aa05aa95SBarry Smith         v = a->v;
623aa05aa95SBarry Smith         for (j=0; j<N; j++) {
624aa05aa95SBarry Smith           for (i=0; i<m; i++) {
625578230a0SSatish Balay             work[j + i*N] = *v++;
626aa05aa95SBarry Smith           }
627aa05aa95SBarry Smith         }
628aa05aa95SBarry Smith         ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
629aa05aa95SBarry Smith         /* get largest work array to receive messages from other processes, excludes process zero */
630aa05aa95SBarry Smith         mmax = 0;
631aa05aa95SBarry Smith         for (i=1; i<size; i++) {
632d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
633aa05aa95SBarry Smith         }
634785e854fSJed Brown         ierr = PetscMalloc1(mmax*N,&vv);CHKERRQ(ierr);
635aa05aa95SBarry Smith         for (k = 1; k < size; k++) {
636f8009846SMatthew Knepley           v    = vv;
637d0f46423SBarry Smith           m    = mat->rmap->range[k+1] - mat->rmap->range[k];
638ce94432eSBarry Smith           ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
639aa05aa95SBarry Smith 
640aa05aa95SBarry Smith           for (j = 0; j < N; j++) {
641aa05aa95SBarry Smith             for (i = 0; i < m; i++) {
642578230a0SSatish Balay               work[j + i*N] = *v++;
643aa05aa95SBarry Smith             }
644aa05aa95SBarry Smith           }
645aa05aa95SBarry Smith           ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
646aa05aa95SBarry Smith         }
647aa05aa95SBarry Smith         ierr = PetscFree(work);CHKERRQ(ierr);
648578230a0SSatish Balay         ierr = PetscFree(vv);CHKERRQ(ierr);
649aa05aa95SBarry Smith       } else {
650ce94432eSBarry Smith         ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
651aa05aa95SBarry Smith       }
652eb3f98d2SBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE)");
653aa05aa95SBarry Smith   }
6543a40ed3dSBarry Smith   PetscFunctionReturn(0);
6558965ea79SLois Curfman McInnes }
6568965ea79SLois Curfman McInnes 
6577da1fb6eSBarry Smith extern PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
6589804daf3SBarry Smith #include <petscdraw.h>
6594a2ae208SSatish Balay #undef __FUNCT__
6604a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
6616849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6628965ea79SLois Curfman McInnes {
66339ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
664dfbe8321SBarry Smith   PetscErrorCode    ierr;
6657da1fb6eSBarry Smith   PetscMPIInt       rank = mdn->rank;
66619fd82e9SBarry Smith   PetscViewerType   vtype;
667ace3abfcSBarry Smith   PetscBool         iascii,isdraw;
668b0a32e0cSBarry Smith   PetscViewer       sviewer;
669f3ef73ceSBarry Smith   PetscViewerFormat format;
6708965ea79SLois Curfman McInnes 
6713a40ed3dSBarry Smith   PetscFunctionBegin;
672251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
673251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
67432077d6dSBarry Smith   if (iascii) {
675b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
676b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
677456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6784e220ebcSLois Curfman McInnes       MatInfo info;
679888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
6801575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
6817b23a99aSBarry 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);
682b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6831575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
6845d00a290SHong Zhang       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
6853a40ed3dSBarry Smith       PetscFunctionReturn(0);
686fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6873a40ed3dSBarry Smith       PetscFunctionReturn(0);
6888965ea79SLois Curfman McInnes     }
689f1af5d2fSBarry Smith   } else if (isdraw) {
690b0a32e0cSBarry Smith     PetscDraw draw;
691ace3abfcSBarry Smith     PetscBool isnull;
692f1af5d2fSBarry Smith 
693b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
694b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
695f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
696f1af5d2fSBarry Smith   }
69777ed5343SBarry Smith 
6987da1fb6eSBarry Smith   {
6998965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
7008965ea79SLois Curfman McInnes     Mat         A;
701d0f46423SBarry Smith     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
702ba8c8a56SBarry Smith     PetscInt    *cols;
703ba8c8a56SBarry Smith     PetscScalar *vals;
7048965ea79SLois Curfman McInnes 
705ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr);
7068965ea79SLois Curfman McInnes     if (!rank) {
707f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
7083a40ed3dSBarry Smith     } else {
709f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
7108965ea79SLois Curfman McInnes     }
7117adad957SLisandro Dalcin     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
712878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
7130298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(A,NULL);CHKERRQ(ierr);
7143bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)A);CHKERRQ(ierr);
7158965ea79SLois Curfman McInnes 
71639ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
71739ddd567SLois Curfman McInnes        but it's quick for now */
71851022da4SBarry Smith     A->insertmode = INSERT_VALUES;
7192205254eSKarl Rupp 
7202205254eSKarl Rupp     row = mat->rmap->rstart;
7212205254eSKarl Rupp     m   = mdn->A->rmap->n;
7228965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
723ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
724ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
725ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
72639ddd567SLois Curfman McInnes       row++;
7278965ea79SLois Curfman McInnes     }
7288965ea79SLois Curfman McInnes 
7296d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7306d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7313f08860eSBarry Smith     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
732b9b97703SBarry Smith     if (!rank) {
7337da1fb6eSBarry Smith         ierr = MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
7348965ea79SLois Curfman McInnes     }
7353f08860eSBarry Smith     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
736b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7376bf464f9SBarry Smith     ierr = MatDestroy(&A);CHKERRQ(ierr);
7388965ea79SLois Curfman McInnes   }
7393a40ed3dSBarry Smith   PetscFunctionReturn(0);
7408965ea79SLois Curfman McInnes }
7418965ea79SLois Curfman McInnes 
7424a2ae208SSatish Balay #undef __FUNCT__
7434a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
744dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
7458965ea79SLois Curfman McInnes {
746dfbe8321SBarry Smith   PetscErrorCode ierr;
747ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw,issocket;
7488965ea79SLois Curfman McInnes 
749433994e6SBarry Smith   PetscFunctionBegin;
750251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
751251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
752251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr);
753251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
7540f5bd95cSBarry Smith 
75532077d6dSBarry Smith   if (iascii || issocket || isdraw) {
756f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7570f5bd95cSBarry Smith   } else if (isbinary) {
7583a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
75911aeaf0aSBarry Smith   }
7603a40ed3dSBarry Smith   PetscFunctionReturn(0);
7618965ea79SLois Curfman McInnes }
7628965ea79SLois Curfman McInnes 
7634a2ae208SSatish Balay #undef __FUNCT__
7644a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
765dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7668965ea79SLois Curfman McInnes {
7673501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7683501a2bdSLois Curfman McInnes   Mat            mdn  = mat->A;
769dfbe8321SBarry Smith   PetscErrorCode ierr;
770329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7718965ea79SLois Curfman McInnes 
7723a40ed3dSBarry Smith   PetscFunctionBegin;
7734e220ebcSLois Curfman McInnes   info->block_size = 1.0;
7742205254eSKarl Rupp 
7754e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7762205254eSKarl Rupp 
7774e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7784e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7798965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7804e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7814e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7824e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7834e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7844e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7858965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
786*b2566f29SBarry Smith     ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
7872205254eSKarl Rupp 
7884e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7894e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7904e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7914e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7924e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7938965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
794*b2566f29SBarry Smith     ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
7952205254eSKarl Rupp 
7964e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7974e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7984e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7994e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
8004e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
8018965ea79SLois Curfman McInnes   }
8024e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
8034e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
8044e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
8053a40ed3dSBarry Smith   PetscFunctionReturn(0);
8068965ea79SLois Curfman McInnes }
8078965ea79SLois Curfman McInnes 
8084a2ae208SSatish Balay #undef __FUNCT__
8094a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
810ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg)
8118965ea79SLois Curfman McInnes {
81239ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
813dfbe8321SBarry Smith   PetscErrorCode ierr;
8148965ea79SLois Curfman McInnes 
8153a40ed3dSBarry Smith   PetscFunctionBegin;
81612c028f9SKris Buschelman   switch (op) {
817512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
81812c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
81912c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
8204e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82112c028f9SKris Buschelman     break;
82212c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
8234e0d8c25SBarry Smith     a->roworiented = flg;
8242205254eSKarl Rupp 
8254e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82612c028f9SKris Buschelman     break;
8274e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
82813fa8e87SLisandro Dalcin   case MAT_KEEP_NONZERO_PATTERN:
82912c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
830290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
83112c028f9SKris Buschelman     break;
83212c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
8334e0d8c25SBarry Smith     a->donotstash = flg;
83412c028f9SKris Buschelman     break;
83577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
83677e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8379a4540c5SBarry Smith   case MAT_HERMITIAN:
8389a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
839600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
840290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
84177e54ba9SKris Buschelman     break;
84212c028f9SKris Buschelman   default:
843e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
8443a40ed3dSBarry Smith   }
8453a40ed3dSBarry Smith   PetscFunctionReturn(0);
8468965ea79SLois Curfman McInnes }
8478965ea79SLois Curfman McInnes 
8488965ea79SLois Curfman McInnes 
8494a2ae208SSatish Balay #undef __FUNCT__
8504a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
851dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8525b2fa520SLois Curfman McInnes {
8535b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8545b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
85587828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
856dfbe8321SBarry Smith   PetscErrorCode ierr;
857d0f46423SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n;
8585b2fa520SLois Curfman McInnes 
8595b2fa520SLois Curfman McInnes   PetscFunctionBegin;
86072d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8615b2fa520SLois Curfman McInnes   if (ll) {
86272d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
863e32f2f54SBarry Smith     if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8641ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8655b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8665b2fa520SLois Curfman McInnes       x = l[i];
8675b2fa520SLois Curfman McInnes       v = mat->v + i;
8685b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8695b2fa520SLois Curfman McInnes     }
8701ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
871efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8725b2fa520SLois Curfman McInnes   }
8735b2fa520SLois Curfman McInnes   if (rr) {
874175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
875e32f2f54SBarry Smith     if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
876ca9f406cSSatish Balay     ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
877ca9f406cSSatish Balay     ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8781ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8795b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8805b2fa520SLois Curfman McInnes       x = r[i];
8815b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8822205254eSKarl Rupp       for (j=0; j<m; j++) (*v++) *= x;
8835b2fa520SLois Curfman McInnes     }
8841ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
885efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8865b2fa520SLois Curfman McInnes   }
8875b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8885b2fa520SLois Curfman McInnes }
8895b2fa520SLois Curfman McInnes 
8904a2ae208SSatish Balay #undef __FUNCT__
8914a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
892dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
893096963f5SLois Curfman McInnes {
8943501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8953501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
896dfbe8321SBarry Smith   PetscErrorCode ierr;
89713f74950SBarry Smith   PetscInt       i,j;
898329f5518SBarry Smith   PetscReal      sum = 0.0;
89987828ca2SBarry Smith   PetscScalar    *v  = mat->v;
9003501a2bdSLois Curfman McInnes 
9013a40ed3dSBarry Smith   PetscFunctionBegin;
9023501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
903064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
9043501a2bdSLois Curfman McInnes   } else {
9053501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
906d0f46423SBarry Smith       for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) {
907329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
9083501a2bdSLois Curfman McInnes       }
909*b2566f29SBarry Smith       ierr = MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
9108f1a2a5eSBarry Smith       *nrm = PetscSqrtReal(*nrm);
911dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr);
9123a40ed3dSBarry Smith     } else if (type == NORM_1) {
913329f5518SBarry Smith       PetscReal *tmp,*tmp2;
914dcca6d9dSJed Brown       ierr = PetscMalloc2(A->cmap->N,&tmp,A->cmap->N,&tmp2);CHKERRQ(ierr);
91574ed9c26SBarry Smith       ierr = PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
91674ed9c26SBarry Smith       ierr = PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
917064f8208SBarry Smith       *nrm = 0.0;
9183501a2bdSLois Curfman McInnes       v    = mat->v;
919d0f46423SBarry Smith       for (j=0; j<mdn->A->cmap->n; j++) {
920d0f46423SBarry Smith         for (i=0; i<mdn->A->rmap->n; i++) {
92167e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
9223501a2bdSLois Curfman McInnes         }
9233501a2bdSLois Curfman McInnes       }
924*b2566f29SBarry Smith       ierr = MPIU_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
925d0f46423SBarry Smith       for (j=0; j<A->cmap->N; j++) {
926064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
9273501a2bdSLois Curfman McInnes       }
9288627564fSBarry Smith       ierr = PetscFree2(tmp,tmp2);CHKERRQ(ierr);
929d0f46423SBarry Smith       ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
9303a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
931329f5518SBarry Smith       PetscReal ntemp;
9323501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
933*b2566f29SBarry Smith       ierr = MPIU_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
934ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm");
9353501a2bdSLois Curfman McInnes   }
9363a40ed3dSBarry Smith   PetscFunctionReturn(0);
9373501a2bdSLois Curfman McInnes }
9383501a2bdSLois Curfman McInnes 
9394a2ae208SSatish Balay #undef __FUNCT__
9404a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
941fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
9423501a2bdSLois Curfman McInnes {
9433501a2bdSLois Curfman McInnes   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
9443501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
9453501a2bdSLois Curfman McInnes   Mat            B;
946d0f46423SBarry Smith   PetscInt       M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart;
9476849ba73SBarry Smith   PetscErrorCode ierr;
94813f74950SBarry Smith   PetscInt       j,i;
94987828ca2SBarry Smith   PetscScalar    *v;
9503501a2bdSLois Curfman McInnes 
9513a40ed3dSBarry Smith   PetscFunctionBegin;
952ce94432eSBarry Smith   if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place");
953fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
954ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
955d0f46423SBarry Smith     ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr);
9567adad957SLisandro Dalcin     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
9570298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr);
958fc4dec0aSBarry Smith   } else {
959fc4dec0aSBarry Smith     B = *matout;
960fc4dec0aSBarry Smith   }
9613501a2bdSLois Curfman McInnes 
962d0f46423SBarry Smith   m    = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v;
963785e854fSJed Brown   ierr = PetscMalloc1(m,&rwork);CHKERRQ(ierr);
9643501a2bdSLois Curfman McInnes   for (i=0; i<m; i++) rwork[i] = rstart + i;
9651acff37aSSatish Balay   for (j=0; j<n; j++) {
9663501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9673501a2bdSLois Curfman McInnes     v   += m;
9683501a2bdSLois Curfman McInnes   }
969606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9706d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9716d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
972815cbec1SBarry Smith   if (reuse == MAT_INITIAL_MATRIX || *matout != A) {
9733501a2bdSLois Curfman McInnes     *matout = B;
9743501a2bdSLois Curfman McInnes   } else {
97528be2f97SBarry Smith     ierr = MatHeaderMerge(A,&B);CHKERRQ(ierr);
9763501a2bdSLois Curfman McInnes   }
9773a40ed3dSBarry Smith   PetscFunctionReturn(0);
978096963f5SLois Curfman McInnes }
979096963f5SLois Curfman McInnes 
98044cd7ae7SLois Curfman McInnes 
9816849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*);
982d84338a6SBarry Smith extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar);
9838965ea79SLois Curfman McInnes 
9844a2ae208SSatish Balay #undef __FUNCT__
9854994cf47SJed Brown #define __FUNCT__ "MatSetUp_MPIDense"
9864994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A)
987273d9f13SBarry Smith {
988dfbe8321SBarry Smith   PetscErrorCode ierr;
989273d9f13SBarry Smith 
990273d9f13SBarry Smith   PetscFunctionBegin;
991273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
992273d9f13SBarry Smith   PetscFunctionReturn(0);
993273d9f13SBarry Smith }
994273d9f13SBarry Smith 
99530716080SHong Zhang #undef __FUNCT__
996488007eeSBarry Smith #define __FUNCT__ "MatAXPY_MPIDense"
997488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
998488007eeSBarry Smith {
999488007eeSBarry Smith   PetscErrorCode ierr;
1000488007eeSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;
1001488007eeSBarry Smith 
1002488007eeSBarry Smith   PetscFunctionBegin;
1003488007eeSBarry Smith   ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr);
1004a3fa217bSJose E. Roman   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
1005488007eeSBarry Smith   PetscFunctionReturn(0);
1006488007eeSBarry Smith }
1007488007eeSBarry Smith 
1008ba337c44SJed Brown #undef __FUNCT__
1009ba337c44SJed Brown #define __FUNCT__ "MatConjugate_MPIDense"
10107087cfbeSBarry Smith PetscErrorCode  MatConjugate_MPIDense(Mat mat)
1011ba337c44SJed Brown {
1012ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)mat->data;
1013ba337c44SJed Brown   PetscErrorCode ierr;
1014ba337c44SJed Brown 
1015ba337c44SJed Brown   PetscFunctionBegin;
1016ba337c44SJed Brown   ierr = MatConjugate(a->A);CHKERRQ(ierr);
1017ba337c44SJed Brown   PetscFunctionReturn(0);
1018ba337c44SJed Brown }
1019ba337c44SJed Brown 
1020ba337c44SJed Brown #undef __FUNCT__
1021ba337c44SJed Brown #define __FUNCT__ "MatRealPart_MPIDense"
1022ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A)
1023ba337c44SJed Brown {
1024ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1025ba337c44SJed Brown   PetscErrorCode ierr;
1026ba337c44SJed Brown 
1027ba337c44SJed Brown   PetscFunctionBegin;
1028ba337c44SJed Brown   ierr = MatRealPart(a->A);CHKERRQ(ierr);
1029ba337c44SJed Brown   PetscFunctionReturn(0);
1030ba337c44SJed Brown }
1031ba337c44SJed Brown 
1032ba337c44SJed Brown #undef __FUNCT__
1033ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_MPIDense"
1034ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1035ba337c44SJed Brown {
1036ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1037ba337c44SJed Brown   PetscErrorCode ierr;
1038ba337c44SJed Brown 
1039ba337c44SJed Brown   PetscFunctionBegin;
1040ba337c44SJed Brown   ierr = MatImaginaryPart(a->A);CHKERRQ(ierr);
1041ba337c44SJed Brown   PetscFunctionReturn(0);
1042ba337c44SJed Brown }
1043ba337c44SJed Brown 
10440716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);
10450716a85fSBarry Smith #undef __FUNCT__
10460716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense"
10470716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
10480716a85fSBarry Smith {
10490716a85fSBarry Smith   PetscErrorCode ierr;
10500716a85fSBarry Smith   PetscInt       i,n;
10510716a85fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
10520716a85fSBarry Smith   PetscReal      *work;
10530716a85fSBarry Smith 
10540716a85fSBarry Smith   PetscFunctionBegin;
10550298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr);
1056785e854fSJed Brown   ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
10570716a85fSBarry Smith   ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr);
10580716a85fSBarry Smith   if (type == NORM_2) {
10590716a85fSBarry Smith     for (i=0; i<n; i++) work[i] *= work[i];
10600716a85fSBarry Smith   }
10610716a85fSBarry Smith   if (type == NORM_INFINITY) {
1062*b2566f29SBarry Smith     ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr);
10630716a85fSBarry Smith   } else {
1064*b2566f29SBarry Smith     ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr);
10650716a85fSBarry Smith   }
10660716a85fSBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
10670716a85fSBarry Smith   if (type == NORM_2) {
10688f1a2a5eSBarry Smith     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
10690716a85fSBarry Smith   }
10700716a85fSBarry Smith   PetscFunctionReturn(0);
10710716a85fSBarry Smith }
10720716a85fSBarry Smith 
107373a71a0fSBarry Smith #undef __FUNCT__
107473a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_MPIDense"
107573a71a0fSBarry Smith static PetscErrorCode  MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
107673a71a0fSBarry Smith {
107773a71a0fSBarry Smith   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;
107873a71a0fSBarry Smith   PetscErrorCode ierr;
107973a71a0fSBarry Smith   PetscScalar    *a;
108073a71a0fSBarry Smith   PetscInt       m,n,i;
108173a71a0fSBarry Smith 
108273a71a0fSBarry Smith   PetscFunctionBegin;
108373a71a0fSBarry Smith   ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr);
10848c778c55SBarry Smith   ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr);
108573a71a0fSBarry Smith   for (i=0; i<m*n; i++) {
108673a71a0fSBarry Smith     ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr);
108773a71a0fSBarry Smith   }
10888c778c55SBarry Smith   ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr);
108973a71a0fSBarry Smith   PetscFunctionReturn(0);
109073a71a0fSBarry Smith }
109173a71a0fSBarry Smith 
1092fd4e9aacSBarry Smith extern PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat,Mat);
1093fd4e9aacSBarry Smith 
10948965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
109509dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense,
109609dc0095SBarry Smith                                         MatGetRow_MPIDense,
109709dc0095SBarry Smith                                         MatRestoreRow_MPIDense,
109809dc0095SBarry Smith                                         MatMult_MPIDense,
109997304618SKris Buschelman                                 /*  4*/ MatMultAdd_MPIDense,
11007c922b88SBarry Smith                                         MatMultTranspose_MPIDense,
11017c922b88SBarry Smith                                         MatMultTransposeAdd_MPIDense,
11028965ea79SLois Curfman McInnes                                         0,
110309dc0095SBarry Smith                                         0,
110409dc0095SBarry Smith                                         0,
110597304618SKris Buschelman                                 /* 10*/ 0,
110609dc0095SBarry Smith                                         0,
110709dc0095SBarry Smith                                         0,
110809dc0095SBarry Smith                                         0,
110909dc0095SBarry Smith                                         MatTranspose_MPIDense,
111097304618SKris Buschelman                                 /* 15*/ MatGetInfo_MPIDense,
11116e4ee0c6SHong Zhang                                         MatEqual_MPIDense,
111209dc0095SBarry Smith                                         MatGetDiagonal_MPIDense,
11135b2fa520SLois Curfman McInnes                                         MatDiagonalScale_MPIDense,
111409dc0095SBarry Smith                                         MatNorm_MPIDense,
111597304618SKris Buschelman                                 /* 20*/ MatAssemblyBegin_MPIDense,
111609dc0095SBarry Smith                                         MatAssemblyEnd_MPIDense,
111709dc0095SBarry Smith                                         MatSetOption_MPIDense,
111809dc0095SBarry Smith                                         MatZeroEntries_MPIDense,
1119d519adbfSMatthew Knepley                                 /* 24*/ MatZeroRows_MPIDense,
1120919b68f7SBarry Smith                                         0,
112101b82886SBarry Smith                                         0,
112201b82886SBarry Smith                                         0,
112301b82886SBarry Smith                                         0,
11244994cf47SJed Brown                                 /* 29*/ MatSetUp_MPIDense,
1125273d9f13SBarry Smith                                         0,
112609dc0095SBarry Smith                                         0,
11278c778c55SBarry Smith                                         0,
11288c778c55SBarry Smith                                         0,
1129d519adbfSMatthew Knepley                                 /* 34*/ MatDuplicate_MPIDense,
113009dc0095SBarry Smith                                         0,
113109dc0095SBarry Smith                                         0,
113209dc0095SBarry Smith                                         0,
113309dc0095SBarry Smith                                         0,
1134d519adbfSMatthew Knepley                                 /* 39*/ MatAXPY_MPIDense,
11352ce60cd0SSatish Balay                                         MatGetSubMatrices_MPIDense,
113609dc0095SBarry Smith                                         0,
113709dc0095SBarry Smith                                         MatGetValues_MPIDense,
113809dc0095SBarry Smith                                         0,
1139d519adbfSMatthew Knepley                                 /* 44*/ 0,
114009dc0095SBarry Smith                                         MatScale_MPIDense,
11417d68702bSBarry Smith                                         MatShift_Basic,
114209dc0095SBarry Smith                                         0,
114309dc0095SBarry Smith                                         0,
114473a71a0fSBarry Smith                                 /* 49*/ MatSetRandom_MPIDense,
114509dc0095SBarry Smith                                         0,
114609dc0095SBarry Smith                                         0,
114709dc0095SBarry Smith                                         0,
114809dc0095SBarry Smith                                         0,
1149d519adbfSMatthew Knepley                                 /* 54*/ 0,
115009dc0095SBarry Smith                                         0,
115109dc0095SBarry Smith                                         0,
115209dc0095SBarry Smith                                         0,
115309dc0095SBarry Smith                                         0,
1154d519adbfSMatthew Knepley                                 /* 59*/ MatGetSubMatrix_MPIDense,
1155b9b97703SBarry Smith                                         MatDestroy_MPIDense,
1156b9b97703SBarry Smith                                         MatView_MPIDense,
1157357abbc8SBarry Smith                                         0,
115897304618SKris Buschelman                                         0,
1159d519adbfSMatthew Knepley                                 /* 64*/ 0,
116097304618SKris Buschelman                                         0,
116197304618SKris Buschelman                                         0,
116297304618SKris Buschelman                                         0,
116397304618SKris Buschelman                                         0,
1164d519adbfSMatthew Knepley                                 /* 69*/ 0,
116597304618SKris Buschelman                                         0,
116697304618SKris Buschelman                                         0,
116797304618SKris Buschelman                                         0,
116897304618SKris Buschelman                                         0,
1169d519adbfSMatthew Knepley                                 /* 74*/ 0,
117097304618SKris Buschelman                                         0,
117197304618SKris Buschelman                                         0,
117297304618SKris Buschelman                                         0,
117397304618SKris Buschelman                                         0,
1174d519adbfSMatthew Knepley                                 /* 79*/ 0,
117597304618SKris Buschelman                                         0,
117697304618SKris Buschelman                                         0,
117797304618SKris Buschelman                                         0,
11785bba2384SShri Abhyankar                                 /* 83*/ MatLoad_MPIDense,
1179865e5f61SKris Buschelman                                         0,
1180865e5f61SKris Buschelman                                         0,
1181865e5f61SKris Buschelman                                         0,
1182865e5f61SKris Buschelman                                         0,
1183865e5f61SKris Buschelman                                         0,
1184d519adbfSMatthew Knepley                                 /* 89*/
1185865e5f61SKris Buschelman                                         0,
1186865e5f61SKris Buschelman                                         0,
1187fd4e9aacSBarry Smith                                         MatMatMultNumeric_MPIDense,
11882fbe02b9SBarry Smith                                         0,
1189ba337c44SJed Brown                                         0,
1190d519adbfSMatthew Knepley                                 /* 94*/ 0,
1191865e5f61SKris Buschelman                                         0,
1192865e5f61SKris Buschelman                                         0,
1193ba337c44SJed Brown                                         0,
1194ba337c44SJed Brown                                         0,
1195ba337c44SJed Brown                                 /* 99*/ 0,
1196ba337c44SJed Brown                                         0,
1197ba337c44SJed Brown                                         0,
1198ba337c44SJed Brown                                         MatConjugate_MPIDense,
1199ba337c44SJed Brown                                         0,
1200ba337c44SJed Brown                                 /*104*/ 0,
1201ba337c44SJed Brown                                         MatRealPart_MPIDense,
1202ba337c44SJed Brown                                         MatImaginaryPart_MPIDense,
120386d161a7SShri Abhyankar                                         0,
120486d161a7SShri Abhyankar                                         0,
120586d161a7SShri Abhyankar                                 /*109*/ 0,
120686d161a7SShri Abhyankar                                         0,
120786d161a7SShri Abhyankar                                         0,
120886d161a7SShri Abhyankar                                         0,
120986d161a7SShri Abhyankar                                         0,
121086d161a7SShri Abhyankar                                 /*114*/ 0,
121186d161a7SShri Abhyankar                                         0,
121286d161a7SShri Abhyankar                                         0,
121386d161a7SShri Abhyankar                                         0,
121486d161a7SShri Abhyankar                                         0,
121586d161a7SShri Abhyankar                                 /*119*/ 0,
121686d161a7SShri Abhyankar                                         0,
121786d161a7SShri Abhyankar                                         0,
12180716a85fSBarry Smith                                         0,
12190716a85fSBarry Smith                                         0,
12200716a85fSBarry Smith                                 /*124*/ 0,
12213964eb88SJed Brown                                         MatGetColumnNorms_MPIDense,
12223964eb88SJed Brown                                         0,
12233964eb88SJed Brown                                         0,
12243964eb88SJed Brown                                         0,
12253964eb88SJed Brown                                 /*129*/ 0,
12263964eb88SJed Brown                                         0,
12273964eb88SJed Brown                                         0,
12283964eb88SJed Brown                                         0,
12293964eb88SJed Brown                                         0,
12303964eb88SJed Brown                                 /*134*/ 0,
12313964eb88SJed Brown                                         0,
12323964eb88SJed Brown                                         0,
12333964eb88SJed Brown                                         0,
12343964eb88SJed Brown                                         0,
12353964eb88SJed Brown                                 /*139*/ 0,
1236f9426fe0SMark Adams                                         0,
12373964eb88SJed Brown                                         0
1238ba337c44SJed Brown };
12398965ea79SLois Curfman McInnes 
12404a2ae208SSatish Balay #undef __FUNCT__
1241a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
12427087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1243a23d5eceSKris Buschelman {
1244a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1245dfbe8321SBarry Smith   PetscErrorCode ierr;
1246a23d5eceSKris Buschelman 
1247a23d5eceSKris Buschelman   PetscFunctionBegin;
1248a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1249a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1250a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1251a23d5eceSKris Buschelman 
1252a23d5eceSKris Buschelman   a       = (Mat_MPIDense*)mat->data;
125334ef9618SShri Abhyankar   ierr    = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
125434ef9618SShri Abhyankar   ierr    = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
125534ef9618SShri Abhyankar   a->nvec = mat->cmap->n;
125634ef9618SShri Abhyankar 
1257f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1258d0f46423SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr);
12595c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
12605c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
12613bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr);
1262a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1263a23d5eceSKris Buschelman }
1264a23d5eceSKris Buschelman 
126565b80a83SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
1266a23d5eceSKris Buschelman #undef __FUNCT__
12678baccfbdSHong Zhang #define __FUNCT__ "MatConvert_MPIDense_Elemental"
12688baccfbdSHong Zhang PETSC_EXTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
12698baccfbdSHong Zhang {
12708ea901baSHong Zhang   Mat            mat_elemental;
12718ea901baSHong Zhang   PetscErrorCode ierr;
127265b80a83SHong Zhang   PetscScalar    *array,*v_rowwise;
12738ea901baSHong Zhang   PetscInt       m=A->rmap->n,N=A->cmap->N,rstart=A->rmap->rstart,i,j,k,*rows,*cols;
12748ea901baSHong Zhang 
12758baccfbdSHong Zhang   PetscFunctionBegin;
127665b80a83SHong Zhang   ierr = PetscMalloc3(m*N,&v_rowwise,m,&rows,N,&cols);CHKERRQ(ierr);
12778ea901baSHong Zhang   ierr = MatDenseGetArray(A,&array);CHKERRQ(ierr);
127865b80a83SHong Zhang   /* convert column-wise array into row-wise v_rowwise, see MatSetValues_Elemental() */
12798ea901baSHong Zhang   k = 0;
12808ea901baSHong Zhang   for (j=0; j<N; j++) {
12818ea901baSHong Zhang     cols[j] = j;
12828ea901baSHong Zhang     for (i=0; i<m; i++) {
128365b80a83SHong Zhang       v_rowwise[i*N+j] = array[k++];
12848ea901baSHong Zhang     }
12858ea901baSHong Zhang   }
12868ea901baSHong Zhang   for (i=0; i<m; i++) {
12878ea901baSHong Zhang     rows[i] = rstart + i;
12888ea901baSHong Zhang   }
12898ea901baSHong Zhang   ierr = MatDenseRestoreArray(A,&array);CHKERRQ(ierr);
12908ea901baSHong Zhang 
12918ea901baSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
12928ea901baSHong Zhang   ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
12938ea901baSHong Zhang   ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
12948ea901baSHong Zhang   ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
12958ea901baSHong Zhang 
12968ea901baSHong Zhang   /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
129765b80a83SHong Zhang   ierr = MatSetValues(mat_elemental,m,rows,N,cols,v_rowwise,ADD_VALUES);CHKERRQ(ierr);
12988ea901baSHong Zhang   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12998ea901baSHong Zhang   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
130065b80a83SHong Zhang   ierr = PetscFree3(v_rowwise,rows,cols);CHKERRQ(ierr);
13018ea901baSHong Zhang 
13028ea901baSHong Zhang   if (reuse == MAT_REUSE_MATRIX) {
130328be2f97SBarry Smith     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
13048ea901baSHong Zhang   } else {
13058ea901baSHong Zhang     *newmat = mat_elemental;
13068ea901baSHong Zhang   }
13078baccfbdSHong Zhang   PetscFunctionReturn(0);
13088baccfbdSHong Zhang }
130965b80a83SHong Zhang #endif
13108baccfbdSHong Zhang 
13118baccfbdSHong Zhang #undef __FUNCT__
13124a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
13138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat)
1314273d9f13SBarry Smith {
1315273d9f13SBarry Smith   Mat_MPIDense   *a;
1316dfbe8321SBarry Smith   PetscErrorCode ierr;
1317273d9f13SBarry Smith 
1318273d9f13SBarry Smith   PetscFunctionBegin;
1319b00a9115SJed Brown   ierr      = PetscNewLog(mat,&a);CHKERRQ(ierr);
1320b0a32e0cSBarry Smith   mat->data = (void*)a;
1321273d9f13SBarry Smith   ierr      = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1322273d9f13SBarry Smith 
1323273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1324ce94432eSBarry Smith   ierr            = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);CHKERRQ(ierr);
1325ce94432eSBarry Smith   ierr            = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);CHKERRQ(ierr);
1326273d9f13SBarry Smith 
1327273d9f13SBarry Smith   /* build cache for off array entries formed */
1328273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
13292205254eSKarl Rupp 
1330ce94432eSBarry Smith   ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);CHKERRQ(ierr);
1331273d9f13SBarry Smith 
1332273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1333273d9f13SBarry Smith   a->lvec        = 0;
1334273d9f13SBarry Smith   a->Mvctx       = 0;
1335273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1336273d9f13SBarry Smith 
1337bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);CHKERRQ(ierr);
1338bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr);
13398baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
13408baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",MatConvert_MPIDense_Elemental);CHKERRQ(ierr);
13418baccfbdSHong Zhang #endif
1342bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1343bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1344bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
1345bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
1346bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
13478949adfdSHong Zhang 
13488949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",MatTransposeMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
13498949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",MatTransposeMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
13508949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",MatTransposeMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
135138aed534SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr);
1352273d9f13SBarry Smith   PetscFunctionReturn(0);
1353273d9f13SBarry Smith }
1354273d9f13SBarry Smith 
1355209238afSKris Buschelman /*MC
1356002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1357209238afSKris Buschelman 
1358209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1359209238afSKris Buschelman    and MATMPIDENSE otherwise.
1360209238afSKris Buschelman 
1361209238afSKris Buschelman    Options Database Keys:
1362209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1363209238afSKris Buschelman 
1364209238afSKris Buschelman   Level: beginner
1365209238afSKris Buschelman 
136601b82886SBarry Smith 
1367209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1368209238afSKris Buschelman M*/
1369209238afSKris Buschelman 
13704a2ae208SSatish Balay #undef __FUNCT__
13714a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1372273d9f13SBarry Smith /*@C
1373273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1374273d9f13SBarry Smith 
1375273d9f13SBarry Smith    Not collective
1376273d9f13SBarry Smith 
1377273d9f13SBarry Smith    Input Parameters:
13781c4f3114SJed Brown .  B - the matrix
13790298fd71SBarry Smith -  data - optional location of matrix data.  Set data=NULL for PETSc
1380273d9f13SBarry Smith    to control all matrix memory allocation.
1381273d9f13SBarry Smith 
1382273d9f13SBarry Smith    Notes:
1383273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1384273d9f13SBarry Smith    storage by columns.
1385273d9f13SBarry Smith 
1386273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1387273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
13880298fd71SBarry Smith    set data=NULL.
1389273d9f13SBarry Smith 
1390273d9f13SBarry Smith    Level: intermediate
1391273d9f13SBarry Smith 
1392273d9f13SBarry Smith .keywords: matrix,dense, parallel
1393273d9f13SBarry Smith 
1394273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1395273d9f13SBarry Smith @*/
13961c4f3114SJed Brown PetscErrorCode  MatMPIDenseSetPreallocation(Mat B,PetscScalar *data)
1397273d9f13SBarry Smith {
13984ac538c5SBarry Smith   PetscErrorCode ierr;
1399273d9f13SBarry Smith 
1400273d9f13SBarry Smith   PetscFunctionBegin;
14011c4f3114SJed Brown   ierr = PetscTryMethod(B,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(B,data));CHKERRQ(ierr);
1402273d9f13SBarry Smith   PetscFunctionReturn(0);
1403273d9f13SBarry Smith }
1404273d9f13SBarry Smith 
14054a2ae208SSatish Balay #undef __FUNCT__
140669b1f4b7SBarry Smith #define __FUNCT__ "MatCreateDense"
14078965ea79SLois Curfman McInnes /*@C
140869b1f4b7SBarry Smith    MatCreateDense - Creates a parallel matrix in dense format.
14098965ea79SLois Curfman McInnes 
1410db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1411db81eaa0SLois Curfman McInnes 
14128965ea79SLois Curfman McInnes    Input Parameters:
1413db81eaa0SLois Curfman McInnes +  comm - MPI communicator
14148965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1415db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
14168965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1417db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
14186cfe35ddSJose E. Roman -  data - optional location of matrix data.  Set data=NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1419dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
14208965ea79SLois Curfman McInnes 
14218965ea79SLois Curfman McInnes    Output Parameter:
1422477f1c0bSLois Curfman McInnes .  A - the matrix
14238965ea79SLois Curfman McInnes 
1424b259b22eSLois Curfman McInnes    Notes:
142539ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
142639ddd567SLois Curfman McInnes    storage by columns.
14278965ea79SLois Curfman McInnes 
142818f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
142918f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
14306cfe35ddSJose E. Roman    set data=NULL (PETSC_NULL_SCALAR for Fortran users).
143118f449edSLois Curfman McInnes 
14328965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
14338965ea79SLois Curfman McInnes    (possibly both).
14348965ea79SLois Curfman McInnes 
1435027ccd11SLois Curfman McInnes    Level: intermediate
1436027ccd11SLois Curfman McInnes 
143739ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
14388965ea79SLois Curfman McInnes 
143939ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
14408965ea79SLois Curfman McInnes @*/
144169b1f4b7SBarry Smith PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
14428965ea79SLois Curfman McInnes {
14436849ba73SBarry Smith   PetscErrorCode ierr;
144413f74950SBarry Smith   PetscMPIInt    size;
14458965ea79SLois Curfman McInnes 
14463a40ed3dSBarry Smith   PetscFunctionBegin;
1447f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1448f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1449273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1450273d9f13SBarry Smith   if (size > 1) {
1451273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1452273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
14536cfe35ddSJose E. Roman     if (data) {  /* user provided data array, so no need to assemble */
14546cfe35ddSJose E. Roman       ierr = MatSetUpMultiply_MPIDense(*A);CHKERRQ(ierr);
14556cfe35ddSJose E. Roman       (*A)->assembled = PETSC_TRUE;
14566cfe35ddSJose E. Roman     }
1457273d9f13SBarry Smith   } else {
1458273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1459273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
14608c469469SLois Curfman McInnes   }
14613a40ed3dSBarry Smith   PetscFunctionReturn(0);
14628965ea79SLois Curfman McInnes }
14638965ea79SLois Curfman McInnes 
14644a2ae208SSatish Balay #undef __FUNCT__
14654a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
14666849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
14678965ea79SLois Curfman McInnes {
14688965ea79SLois Curfman McInnes   Mat            mat;
14693501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1470dfbe8321SBarry Smith   PetscErrorCode ierr;
14718965ea79SLois Curfman McInnes 
14723a40ed3dSBarry Smith   PetscFunctionBegin;
14738965ea79SLois Curfman McInnes   *newmat = 0;
1474ce94432eSBarry Smith   ierr    = MatCreate(PetscObjectComm((PetscObject)A),&mat);CHKERRQ(ierr);
1475d0f46423SBarry Smith   ierr    = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
14767adad957SLisandro Dalcin   ierr    = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1477834f8fabSBarry Smith   a       = (Mat_MPIDense*)mat->data;
1478e04c1aa4SHong Zhang   ierr    = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
14795aa7edbeSHong Zhang 
1480d5f3da31SBarry Smith   mat->factortype   = A->factortype;
1481c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1482273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
14838965ea79SLois Curfman McInnes 
14848965ea79SLois Curfman McInnes   a->size         = oldmat->size;
14858965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1486e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1487b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
14883782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1489e04c1aa4SHong Zhang 
14901e1e43feSBarry Smith   ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr);
14911e1e43feSBarry Smith   ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr);
14928965ea79SLois Curfman McInnes 
1493329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
14945609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
14953bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr);
149601b82886SBarry Smith 
14978965ea79SLois Curfman McInnes   *newmat = mat;
14983a40ed3dSBarry Smith   PetscFunctionReturn(0);
14998965ea79SLois Curfman McInnes }
15008965ea79SLois Curfman McInnes 
15014a2ae208SSatish Balay #undef __FUNCT__
15025bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
15039d36ed5fSBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat)
150486d161a7SShri Abhyankar {
150586d161a7SShri Abhyankar   PetscErrorCode ierr;
150686d161a7SShri Abhyankar   PetscMPIInt    rank,size;
150774c13d95SBarry Smith   const PetscInt *rowners;
150874c13d95SBarry Smith   PetscInt       i,m,n,nz,j,mMax;
150986d161a7SShri Abhyankar   PetscScalar    *array,*vals,*vals_ptr;
15109d36ed5fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)newmat->data;
151186d161a7SShri Abhyankar 
151286d161a7SShri Abhyankar   PetscFunctionBegin;
151386d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
151486d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
151586d161a7SShri Abhyankar 
151674c13d95SBarry Smith   /* determine ownership of rows and columns */
151774c13d95SBarry Smith   m = (newmat->rmap->n < 0) ? PETSC_DECIDE : newmat->rmap->n;
151874c13d95SBarry Smith   n = (newmat->cmap->n < 0) ? PETSC_DECIDE : newmat->cmap->n;
151986d161a7SShri Abhyankar 
152074c13d95SBarry Smith   ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr);
15219d36ed5fSBarry Smith   if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) {
15220298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
15239d36ed5fSBarry Smith   }
15248c778c55SBarry Smith   ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr);
152549467e7dSBarry Smith   ierr = MatGetLocalSize(newmat,&m,NULL);CHKERRQ(ierr);
152674c13d95SBarry Smith   ierr = MatGetOwnershipRanges(newmat,&rowners);CHKERRQ(ierr);
1527396e826eSBarry Smith   ierr = MPI_Reduce(&m,&mMax,1,MPIU_INT,MPI_MAX,0,comm);CHKERRQ(ierr);
152886d161a7SShri Abhyankar   if (!rank) {
15299d36ed5fSBarry Smith     ierr = PetscMalloc1(mMax*N,&vals);CHKERRQ(ierr);
153086d161a7SShri Abhyankar 
153186d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
153286d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
153386d161a7SShri Abhyankar 
153486d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
153586d161a7SShri Abhyankar     vals_ptr = vals;
153686d161a7SShri Abhyankar     for (i=0; i<m; i++) {
153786d161a7SShri Abhyankar       for (j=0; j<N; j++) {
153886d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
153986d161a7SShri Abhyankar       }
154086d161a7SShri Abhyankar     }
154186d161a7SShri Abhyankar 
154286d161a7SShri Abhyankar     /* read in other processors and ship out */
154386d161a7SShri Abhyankar     for (i=1; i<size; i++) {
154486d161a7SShri Abhyankar       nz   = (rowners[i+1] - rowners[i])*N;
154586d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1546a25532f0SBarry Smith       ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
154786d161a7SShri Abhyankar     }
154886d161a7SShri Abhyankar   } else {
154986d161a7SShri Abhyankar     /* receive numeric values */
1550785e854fSJed Brown     ierr = PetscMalloc1(m*N,&vals);CHKERRQ(ierr);
155186d161a7SShri Abhyankar 
155286d161a7SShri Abhyankar     /* receive message of values*/
1553a25532f0SBarry Smith     ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
155486d161a7SShri Abhyankar 
155586d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
155686d161a7SShri Abhyankar     vals_ptr = vals;
155786d161a7SShri Abhyankar     for (i=0; i<m; i++) {
155886d161a7SShri Abhyankar       for (j=0; j<N; j++) {
155986d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
156086d161a7SShri Abhyankar       }
156186d161a7SShri Abhyankar     }
156286d161a7SShri Abhyankar   }
15638c778c55SBarry Smith   ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr);
156486d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
156586d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156686d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156786d161a7SShri Abhyankar   PetscFunctionReturn(0);
156886d161a7SShri Abhyankar }
156986d161a7SShri Abhyankar 
157086d161a7SShri Abhyankar #undef __FUNCT__
15715bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense"
1572112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer)
157386d161a7SShri Abhyankar {
1574dfdea288SBarry Smith   Mat_MPIDense   *a;
157586d161a7SShri Abhyankar   PetscScalar    *vals,*svals;
1576ce94432eSBarry Smith   MPI_Comm       comm;
157786d161a7SShri Abhyankar   MPI_Status     status;
157849467e7dSBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,n,maxnz;
157986d161a7SShri Abhyankar   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
158049467e7dSBarry Smith   PetscInt       *ourlens,*procsnz = 0,jj,*mycols,*smycols;
15819d36ed5fSBarry Smith   PetscInt       i,nz,j,rstart,rend;
158286d161a7SShri Abhyankar   int            fd;
158386d161a7SShri Abhyankar   PetscErrorCode ierr;
158486d161a7SShri Abhyankar 
158586d161a7SShri Abhyankar   PetscFunctionBegin;
1586c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
1587c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
1588ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
158986d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
159086d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
159186d161a7SShri Abhyankar   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
15925872f025SBarry Smith   if (!rank) {
159386d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,(char*)header,4,PETSC_INT);CHKERRQ(ierr);
159486d161a7SShri Abhyankar     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
159586d161a7SShri Abhyankar   }
159686d161a7SShri Abhyankar   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
159786d161a7SShri Abhyankar   M    = header[1]; N = header[2]; nz = header[3];
159886d161a7SShri Abhyankar 
159986d161a7SShri Abhyankar   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
16009d36ed5fSBarry Smith   if (newmat->rmap->N < 0) newmat->rmap->N = M;
16019d36ed5fSBarry Smith   if (newmat->cmap->N < 0) newmat->cmap->N = N;
160286d161a7SShri Abhyankar 
16039d36ed5fSBarry 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);
16049d36ed5fSBarry 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);
160586d161a7SShri Abhyankar 
160686d161a7SShri Abhyankar   /*
160786d161a7SShri Abhyankar        Handle case where matrix is stored on disk as a dense matrix
160886d161a7SShri Abhyankar   */
160986d161a7SShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
16109d36ed5fSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
161186d161a7SShri Abhyankar     PetscFunctionReturn(0);
161286d161a7SShri Abhyankar   }
161386d161a7SShri Abhyankar 
161486d161a7SShri Abhyankar   /* determine ownership of all rows */
16152205254eSKarl Rupp   if (newmat->rmap->n < 0) {
16162205254eSKarl Rupp     ierr = PetscMPIIntCast(M/size + ((M % size) > rank),&m);CHKERRQ(ierr);
16172205254eSKarl Rupp   } else {
16182205254eSKarl Rupp     ierr = PetscMPIIntCast(newmat->rmap->n,&m);CHKERRQ(ierr);
16192205254eSKarl Rupp   }
162049467e7dSBarry Smith   if (newmat->cmap->n < 0) {
162149467e7dSBarry Smith     n = PETSC_DECIDE;
162249467e7dSBarry Smith   } else {
162349467e7dSBarry Smith     ierr = PetscMPIIntCast(newmat->cmap->n,&n);CHKERRQ(ierr);
162449467e7dSBarry Smith   }
162549467e7dSBarry Smith 
1626854ce69bSBarry Smith   ierr       = PetscMalloc1(size+2,&rowners);CHKERRQ(ierr);
162786d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
162886d161a7SShri Abhyankar   rowners[0] = 0;
162986d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
163086d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
163186d161a7SShri Abhyankar   }
163286d161a7SShri Abhyankar   rstart = rowners[rank];
163386d161a7SShri Abhyankar   rend   = rowners[rank+1];
163486d161a7SShri Abhyankar 
163586d161a7SShri Abhyankar   /* distribute row lengths to all processors */
163649467e7dSBarry Smith   ierr = PetscMalloc1(rend-rstart,&ourlens);CHKERRQ(ierr);
163786d161a7SShri Abhyankar   if (!rank) {
1638785e854fSJed Brown     ierr = PetscMalloc1(M,&rowlengths);CHKERRQ(ierr);
163986d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1640785e854fSJed Brown     ierr = PetscMalloc1(size,&sndcounts);CHKERRQ(ierr);
164186d161a7SShri Abhyankar     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
164286d161a7SShri Abhyankar     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
164386d161a7SShri Abhyankar     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
164486d161a7SShri Abhyankar   } else {
164586d161a7SShri Abhyankar     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
164686d161a7SShri Abhyankar   }
164786d161a7SShri Abhyankar 
164886d161a7SShri Abhyankar   if (!rank) {
164986d161a7SShri Abhyankar     /* calculate the number of nonzeros on each processor */
1650785e854fSJed Brown     ierr = PetscMalloc1(size,&procsnz);CHKERRQ(ierr);
165186d161a7SShri Abhyankar     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
165286d161a7SShri Abhyankar     for (i=0; i<size; i++) {
165386d161a7SShri Abhyankar       for (j=rowners[i]; j< rowners[i+1]; j++) {
165486d161a7SShri Abhyankar         procsnz[i] += rowlengths[j];
165586d161a7SShri Abhyankar       }
165686d161a7SShri Abhyankar     }
165786d161a7SShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
165886d161a7SShri Abhyankar 
165986d161a7SShri Abhyankar     /* determine max buffer needed and allocate it */
166086d161a7SShri Abhyankar     maxnz = 0;
166186d161a7SShri Abhyankar     for (i=0; i<size; i++) {
166286d161a7SShri Abhyankar       maxnz = PetscMax(maxnz,procsnz[i]);
166386d161a7SShri Abhyankar     }
1664785e854fSJed Brown     ierr = PetscMalloc1(maxnz,&cols);CHKERRQ(ierr);
166586d161a7SShri Abhyankar 
166686d161a7SShri Abhyankar     /* read in my part of the matrix column indices  */
166786d161a7SShri Abhyankar     nz   = procsnz[0];
1668785e854fSJed Brown     ierr = PetscMalloc1(nz,&mycols);CHKERRQ(ierr);
166986d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
167086d161a7SShri Abhyankar 
167186d161a7SShri Abhyankar     /* read in every one elses and ship off */
167286d161a7SShri Abhyankar     for (i=1; i<size; i++) {
167386d161a7SShri Abhyankar       nz   = procsnz[i];
167486d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
167586d161a7SShri Abhyankar       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
167686d161a7SShri Abhyankar     }
167786d161a7SShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
167886d161a7SShri Abhyankar   } else {
167986d161a7SShri Abhyankar     /* determine buffer space needed for message */
168086d161a7SShri Abhyankar     nz = 0;
168186d161a7SShri Abhyankar     for (i=0; i<m; i++) {
168286d161a7SShri Abhyankar       nz += ourlens[i];
168386d161a7SShri Abhyankar     }
1684854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&mycols);CHKERRQ(ierr);
168586d161a7SShri Abhyankar 
168686d161a7SShri Abhyankar     /* receive message of column indices*/
168786d161a7SShri Abhyankar     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
168886d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
168986d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
169086d161a7SShri Abhyankar   }
169186d161a7SShri Abhyankar 
169249467e7dSBarry Smith   ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr);
1693dfdea288SBarry Smith   a = (Mat_MPIDense*)newmat->data;
16942e3566b0SBarry Smith   if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) {
16950298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
1696dfdea288SBarry Smith   }
169786d161a7SShri Abhyankar 
169886d161a7SShri Abhyankar   if (!rank) {
1699785e854fSJed Brown     ierr = PetscMalloc1(maxnz,&vals);CHKERRQ(ierr);
170086d161a7SShri Abhyankar 
170186d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
170286d161a7SShri Abhyankar     nz   = procsnz[0];
170386d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
170486d161a7SShri Abhyankar 
170586d161a7SShri Abhyankar     /* insert into matrix */
170686d161a7SShri Abhyankar     jj      = rstart;
170786d161a7SShri Abhyankar     smycols = mycols;
170886d161a7SShri Abhyankar     svals   = vals;
170986d161a7SShri Abhyankar     for (i=0; i<m; i++) {
171086d161a7SShri Abhyankar       ierr     = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
171186d161a7SShri Abhyankar       smycols += ourlens[i];
171286d161a7SShri Abhyankar       svals   += ourlens[i];
171386d161a7SShri Abhyankar       jj++;
171486d161a7SShri Abhyankar     }
171586d161a7SShri Abhyankar 
171686d161a7SShri Abhyankar     /* read in other processors and ship out */
171786d161a7SShri Abhyankar     for (i=1; i<size; i++) {
171886d161a7SShri Abhyankar       nz   = procsnz[i];
171986d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
172086d161a7SShri Abhyankar       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr);
172186d161a7SShri Abhyankar     }
172286d161a7SShri Abhyankar     ierr = PetscFree(procsnz);CHKERRQ(ierr);
172386d161a7SShri Abhyankar   } else {
172486d161a7SShri Abhyankar     /* receive numeric values */
1725854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&vals);CHKERRQ(ierr);
172686d161a7SShri Abhyankar 
172786d161a7SShri Abhyankar     /* receive message of values*/
172886d161a7SShri Abhyankar     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr);
172986d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
173086d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
173186d161a7SShri Abhyankar 
173286d161a7SShri Abhyankar     /* insert into matrix */
173386d161a7SShri Abhyankar     jj      = rstart;
173486d161a7SShri Abhyankar     smycols = mycols;
173586d161a7SShri Abhyankar     svals   = vals;
173686d161a7SShri Abhyankar     for (i=0; i<m; i++) {
173786d161a7SShri Abhyankar       ierr     = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
173886d161a7SShri Abhyankar       smycols += ourlens[i];
173986d161a7SShri Abhyankar       svals   += ourlens[i];
174086d161a7SShri Abhyankar       jj++;
174186d161a7SShri Abhyankar     }
174286d161a7SShri Abhyankar   }
174349467e7dSBarry Smith   ierr = PetscFree(ourlens);CHKERRQ(ierr);
174486d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
174586d161a7SShri Abhyankar   ierr = PetscFree(mycols);CHKERRQ(ierr);
174686d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
174786d161a7SShri Abhyankar 
174886d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
174986d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
175086d161a7SShri Abhyankar   PetscFunctionReturn(0);
175186d161a7SShri Abhyankar }
175286d161a7SShri Abhyankar 
175386d161a7SShri Abhyankar #undef __FUNCT__
17546e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense"
1755ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
17566e4ee0c6SHong Zhang {
17576e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
17586e4ee0c6SHong Zhang   Mat            a,b;
1759ace3abfcSBarry Smith   PetscBool      flg;
17606e4ee0c6SHong Zhang   PetscErrorCode ierr;
176190ace30eSBarry Smith 
17626e4ee0c6SHong Zhang   PetscFunctionBegin;
17636e4ee0c6SHong Zhang   a    = matA->A;
17646e4ee0c6SHong Zhang   b    = matB->A;
17656e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
1766*b2566f29SBarry Smith   ierr = MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
17676e4ee0c6SHong Zhang   PetscFunctionReturn(0);
17686e4ee0c6SHong Zhang }
176990ace30eSBarry Smith 
1770