xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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*/
88965ea79SLois Curfman McInnes 
9ba8c8a56SBarry Smith #undef __FUNCT__
10ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
11ab92ecdeSBarry Smith /*@
12ab92ecdeSBarry Smith 
13ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
14ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
15ab92ecdeSBarry Smith 
16ab92ecdeSBarry Smith     Input Parameter:
17ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
18ab92ecdeSBarry Smith 
19ab92ecdeSBarry Smith     Output Parameter:
20ab92ecdeSBarry Smith .      B - the inner matrix
21ab92ecdeSBarry Smith 
228e6c10adSSatish Balay     Level: intermediate
238e6c10adSSatish Balay 
24ab92ecdeSBarry Smith @*/
25ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
26ab92ecdeSBarry Smith {
27ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
28ab92ecdeSBarry Smith   PetscErrorCode ierr;
29ace3abfcSBarry Smith   PetscBool      flg;
30ab92ecdeSBarry Smith 
31ab92ecdeSBarry Smith   PetscFunctionBegin;
32251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
332205254eSKarl Rupp   if (flg) *B = mat->A;
342205254eSKarl Rupp   else *B = A;
35ab92ecdeSBarry Smith   PetscFunctionReturn(0);
36ab92ecdeSBarry Smith }
37ab92ecdeSBarry Smith 
38ab92ecdeSBarry Smith #undef __FUNCT__
39ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
40ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
41ba8c8a56SBarry Smith {
42ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
43ba8c8a56SBarry Smith   PetscErrorCode ierr;
44d0f46423SBarry Smith   PetscInt       lrow,rstart = A->rmap->rstart,rend = A->rmap->rend;
45ba8c8a56SBarry Smith 
46ba8c8a56SBarry Smith   PetscFunctionBegin;
47e7e72b3dSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows");
48ba8c8a56SBarry Smith   lrow = row - rstart;
49ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt**)idx,(const PetscScalar**)v);CHKERRQ(ierr);
50ba8c8a56SBarry Smith   PetscFunctionReturn(0);
51ba8c8a56SBarry Smith }
52ba8c8a56SBarry Smith 
53ba8c8a56SBarry Smith #undef __FUNCT__
54ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
55ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
56ba8c8a56SBarry Smith {
57ba8c8a56SBarry Smith   PetscErrorCode ierr;
58ba8c8a56SBarry Smith 
59ba8c8a56SBarry Smith   PetscFunctionBegin;
60ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
61ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
62ba8c8a56SBarry Smith   PetscFunctionReturn(0);
63ba8c8a56SBarry Smith }
64ba8c8a56SBarry Smith 
650de54da6SSatish Balay EXTERN_C_BEGIN
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 EXTERN_C_END
980de54da6SSatish Balay 
994a2ae208SSatish Balay #undef __FUNCT__
1004a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
10113f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1028965ea79SLois Curfman McInnes {
10339b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
104dfbe8321SBarry Smith   PetscErrorCode ierr;
105d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
106ace3abfcSBarry Smith   PetscBool      roworiented = A->roworiented;
1078965ea79SLois Curfman McInnes 
1083a40ed3dSBarry Smith   PetscFunctionBegin;
10971fd2e92SBarry Smith   if (v) PetscValidScalarPointer(v,6);
1108965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1115ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
112e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1138965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1148965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
11539b7565bSBarry Smith       if (roworiented) {
11639b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1173a40ed3dSBarry Smith       } else {
1188965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1195ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
120e32f2f54SBarry Smith           if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
12139b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
12239b7565bSBarry Smith         }
1238965ea79SLois Curfman McInnes       }
1242205254eSKarl Rupp     } else if (!A->donotstash) {
1255080c13bSMatthew G Knepley       mat->assembled = PETSC_FALSE;
12639b7565bSBarry Smith       if (roworiented) {
127b400d20cSBarry Smith         ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n,PETSC_FALSE);CHKERRQ(ierr);
128d36fbae8SSatish Balay       } else {
129b400d20cSBarry Smith         ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m,PETSC_FALSE);CHKERRQ(ierr);
13039b7565bSBarry Smith       }
131b49de8d1SLois Curfman McInnes     }
132b49de8d1SLois Curfman McInnes   }
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
134b49de8d1SLois Curfman McInnes }
135b49de8d1SLois Curfman McInnes 
1364a2ae208SSatish Balay #undef __FUNCT__
1374a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
13813f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
139b49de8d1SLois Curfman McInnes {
140b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
141dfbe8321SBarry Smith   PetscErrorCode ierr;
142d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
143b49de8d1SLois Curfman McInnes 
1443a40ed3dSBarry Smith   PetscFunctionBegin;
145b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
146e32f2f54SBarry Smith     if (idxm[i] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
147e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
148b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
149b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
150b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
151e32f2f54SBarry Smith         if (idxn[j] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
152e7e72b3dSBarry Smith         if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
153b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
154b49de8d1SLois Curfman McInnes       }
155e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
1568965ea79SLois Curfman McInnes   }
1573a40ed3dSBarry Smith   PetscFunctionReturn(0);
1588965ea79SLois Curfman McInnes }
1598965ea79SLois Curfman McInnes 
1604a2ae208SSatish Balay #undef __FUNCT__
1618c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_MPIDense"
1628c778c55SBarry Smith PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[])
163ff14e315SSatish Balay {
164ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
165dfbe8321SBarry Smith   PetscErrorCode ierr;
166ff14e315SSatish Balay 
1673a40ed3dSBarry Smith   PetscFunctionBegin;
1688c778c55SBarry Smith   ierr = MatDenseGetArray(a->A,array);CHKERRQ(ierr);
1693a40ed3dSBarry Smith   PetscFunctionReturn(0);
170ff14e315SSatish Balay }
171ff14e315SSatish Balay 
1724a2ae208SSatish Balay #undef __FUNCT__
1734a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
1744aa3045dSJed Brown static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
175ca3fa75bSLois Curfman McInnes {
176ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat  = (Mat_MPIDense*)A->data,*newmatd;
177ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
1786849ba73SBarry Smith   PetscErrorCode ierr;
1794aa3045dSJed Brown   PetscInt       i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols;
1805d0c19d7SBarry Smith   const PetscInt *irow,*icol;
18187828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
182ca3fa75bSLois Curfman McInnes   Mat            newmat;
1834aa3045dSJed Brown   IS             iscol_local;
184ca3fa75bSLois Curfman McInnes 
185ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
1864aa3045dSJed Brown   ierr = ISAllGather(iscol,&iscol_local);CHKERRQ(ierr);
187ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
1884aa3045dSJed Brown   ierr = ISGetIndices(iscol_local,&icol);CHKERRQ(ierr);
189b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
190b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
1914aa3045dSJed Brown   ierr = ISGetSize(iscol,&Ncols);CHKERRQ(ierr); /* global number of columns, size of iscol_local */
192ca3fa75bSLois Curfman McInnes 
193ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1947eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1957eba5e9cSLois Curfman McInnes      original matrix! */
196ca3fa75bSLois Curfman McInnes 
197ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1987eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
199ca3fa75bSLois Curfman McInnes 
200ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
201ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
202e32f2f54SBarry Smith     /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2037eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
204ca3fa75bSLois Curfman McInnes     newmat = *B;
205ca3fa75bSLois Curfman McInnes   } else {
206ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
207*ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr);
2084aa3045dSJed Brown     ierr = MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);CHKERRQ(ierr);
2097adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
2100298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
211ca3fa75bSLois Curfman McInnes   }
212ca3fa75bSLois Curfman McInnes 
213ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
214ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
215ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense*)newmatd->A->data)->v;
216ca3fa75bSLois Curfman McInnes 
2174aa3045dSJed Brown   for (i=0; i<Ncols; i++) {
21825a33276SHong Zhang     av = v + ((Mat_SeqDense*)mat->A->data)->lda*icol[i];
219ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2207eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
221ca3fa75bSLois Curfman McInnes     }
222ca3fa75bSLois Curfman McInnes   }
223ca3fa75bSLois Curfman McInnes 
224ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
225ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
226ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
227ca3fa75bSLois Curfman McInnes 
228ca3fa75bSLois Curfman McInnes   /* Free work space */
229ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
2305bdf786aSShri Abhyankar   ierr = ISRestoreIndices(iscol_local,&icol);CHKERRQ(ierr);
23132bb1f2dSLisandro Dalcin   ierr = ISDestroy(&iscol_local);CHKERRQ(ierr);
232ca3fa75bSLois Curfman McInnes   *B   = newmat;
233ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
234ca3fa75bSLois Curfman McInnes }
235ca3fa75bSLois Curfman McInnes 
2364a2ae208SSatish Balay #undef __FUNCT__
2378c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_MPIDense"
2388c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[])
239ff14e315SSatish Balay {
24073a71a0fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
24173a71a0fSBarry Smith   PetscErrorCode ierr;
24273a71a0fSBarry Smith 
2433a40ed3dSBarry Smith   PetscFunctionBegin;
2448c778c55SBarry Smith   ierr = MatDenseRestoreArray(a->A,array);CHKERRQ(ierr);
2453a40ed3dSBarry Smith   PetscFunctionReturn(0);
246ff14e315SSatish Balay }
247ff14e315SSatish Balay 
2484a2ae208SSatish Balay #undef __FUNCT__
2494a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
250dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2518965ea79SLois Curfman McInnes {
25239ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
253*ce94432eSBarry Smith   MPI_Comm       comm;
254dfbe8321SBarry Smith   PetscErrorCode ierr;
25513f74950SBarry Smith   PetscInt       nstash,reallocs;
2568965ea79SLois Curfman McInnes   InsertMode     addv;
2578965ea79SLois Curfman McInnes 
2583a40ed3dSBarry Smith   PetscFunctionBegin;
259*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
2608965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
261c3aae356SJed Brown   ierr = MPI_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);CHKERRQ(ierr);
262*ce94432eSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
263e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2648965ea79SLois Curfman McInnes 
265d0f46423SBarry Smith   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr);
2668798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
267ae15b995SBarry Smith   ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
2683a40ed3dSBarry Smith   PetscFunctionReturn(0);
2698965ea79SLois Curfman McInnes }
2708965ea79SLois Curfman McInnes 
2714a2ae208SSatish Balay #undef __FUNCT__
2724a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
273dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2748965ea79SLois Curfman McInnes {
27539ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn=(Mat_MPIDense*)mat->data;
2766849ba73SBarry Smith   PetscErrorCode ierr;
27713f74950SBarry Smith   PetscInt       i,*row,*col,flg,j,rstart,ncols;
27813f74950SBarry Smith   PetscMPIInt    n;
27987828ca2SBarry Smith   PetscScalar    *val;
280e0fa3b82SLois Curfman McInnes   InsertMode     addv=mat->insertmode;
2818965ea79SLois Curfman McInnes 
2823a40ed3dSBarry Smith   PetscFunctionBegin;
2838965ea79SLois Curfman McInnes   /*  wait on receives */
2847ef1d9bdSSatish Balay   while (1) {
2858798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2867ef1d9bdSSatish Balay     if (!flg) break;
2878965ea79SLois Curfman McInnes 
2887ef1d9bdSSatish Balay     for (i=0; i<n;) {
2897ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2902205254eSKarl Rupp       for (j=i,rstart=row[j]; j<n; j++) {
2912205254eSKarl Rupp         if (row[j] != rstart) break;
2922205254eSKarl Rupp       }
2937ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2947ef1d9bdSSatish Balay       else       ncols = n-i;
2957ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2967ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2977ef1d9bdSSatish Balay       i    = j;
2988965ea79SLois Curfman McInnes     }
2997ef1d9bdSSatish Balay   }
3008798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
3018965ea79SLois Curfman McInnes 
30239ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
30339ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3048965ea79SLois Curfman McInnes 
3056d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
30639ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3078965ea79SLois Curfman McInnes   }
3083a40ed3dSBarry Smith   PetscFunctionReturn(0);
3098965ea79SLois Curfman McInnes }
3108965ea79SLois Curfman McInnes 
3114a2ae208SSatish Balay #undef __FUNCT__
3124a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
313dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3148965ea79SLois Curfman McInnes {
315dfbe8321SBarry Smith   PetscErrorCode ierr;
31639ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3173a40ed3dSBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
3193a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3203a40ed3dSBarry Smith   PetscFunctionReturn(0);
3218965ea79SLois Curfman McInnes }
3228965ea79SLois Curfman McInnes 
3238965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3248965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3258965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3263501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3278965ea79SLois Curfman McInnes    routine.
3288965ea79SLois Curfman McInnes */
3294a2ae208SSatish Balay #undef __FUNCT__
3304a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
3312b40b63fSBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
3328965ea79SLois Curfman McInnes {
33339ddd567SLois Curfman McInnes   Mat_MPIDense      *l = (Mat_MPIDense*)A->data;
3346849ba73SBarry Smith   PetscErrorCode    ierr;
335d0f46423SBarry Smith   PetscInt          i,*owners = A->rmap->range;
33613f74950SBarry Smith   PetscInt          *nprocs,j,idx,nsends;
33713f74950SBarry Smith   PetscInt          nmax,*svalues,*starts,*owner,nrecvs;
3387adad957SLisandro Dalcin   PetscInt          *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source;
33913f74950SBarry Smith   PetscInt          *lens,*lrows,*values;
34013f74950SBarry Smith   PetscMPIInt       n,imdex,rank = l->rank,size = l->size;
341*ce94432eSBarry Smith   MPI_Comm          comm;
3428965ea79SLois Curfman McInnes   MPI_Request       *send_waits,*recv_waits;
3438965ea79SLois Curfman McInnes   MPI_Status        recv_status,*send_status;
344ace3abfcSBarry Smith   PetscBool         found;
34597b48c8fSBarry Smith   const PetscScalar *xx;
34697b48c8fSBarry Smith   PetscScalar       *bb;
3478965ea79SLois Curfman McInnes 
3483a40ed3dSBarry Smith   PetscFunctionBegin;
349*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
350*ce94432eSBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(comm,PETSC_ERR_SUP,"Only handles square matrices");
351b9679d65SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership");
3528965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
35313f74950SBarry Smith   ierr = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
35413f74950SBarry Smith   ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
35513f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr);  /* see note*/
3568965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3578965ea79SLois Curfman McInnes     idx   = rows[i];
35835d8aa7fSBarry Smith     found = PETSC_FALSE;
3598965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3608965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
361c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3628965ea79SLois Curfman McInnes       }
3638965ea79SLois Curfman McInnes     }
364e32f2f54SBarry Smith     if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3658965ea79SLois Curfman McInnes   }
3662205254eSKarl Rupp   nsends = 0;
3672205254eSKarl Rupp   for (i=0; i<size; i++) nsends += nprocs[2*i+1];
3688965ea79SLois Curfman McInnes 
3698965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
370c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3718965ea79SLois Curfman McInnes 
3728965ea79SLois Curfman McInnes   /* post receives:   */
37313f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
374b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3758965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
37613f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes   }
3788965ea79SLois Curfman McInnes 
3798965ea79SLois Curfman McInnes   /* do sends:
3808965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3818965ea79SLois Curfman McInnes          the ith processor
3828965ea79SLois Curfman McInnes   */
38313f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
384b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
38513f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
3868965ea79SLois Curfman McInnes 
3878965ea79SLois Curfman McInnes   starts[0] = 0;
3882205254eSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + nprocs[2*i-2];
3892205254eSKarl Rupp   for (i=0; i<N; i++) svalues[starts[owner[i]]++] = rows[i];
3902205254eSKarl Rupp 
3912205254eSKarl Rupp   starts[0] = 0;
3922205254eSKarl Rupp   for (i=1; i<size+1; i++) starts[i] = starts[i-1] + nprocs[2*i-2];
3938965ea79SLois Curfman McInnes   count = 0;
3948965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
395c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
39613f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3978965ea79SLois Curfman McInnes     }
3988965ea79SLois Curfman McInnes   }
399606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
4008965ea79SLois Curfman McInnes 
4018965ea79SLois Curfman McInnes   base = owners[rank];
4028965ea79SLois Curfman McInnes 
4038965ea79SLois Curfman McInnes   /*  wait on receives */
40474ed9c26SBarry Smith   ierr  = PetscMalloc2(nrecvs,PetscInt,&lens,nrecvs,PetscInt,&source);CHKERRQ(ierr);
40574ed9c26SBarry Smith   count = nrecvs;
40674ed9c26SBarry Smith   slen  = 0;
4078965ea79SLois Curfman McInnes   while (count) {
408ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4098965ea79SLois Curfman McInnes     /* unpack receives into our local space */
41013f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4112205254eSKarl Rupp 
4128965ea79SLois Curfman McInnes     source[imdex] = recv_status.MPI_SOURCE;
4138965ea79SLois Curfman McInnes     lens[imdex]   = n;
4148965ea79SLois Curfman McInnes     slen += n;
4158965ea79SLois Curfman McInnes     count--;
4168965ea79SLois Curfman McInnes   }
417606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4188965ea79SLois Curfman McInnes 
4198965ea79SLois Curfman McInnes   /* move the data into the send scatter */
42013f74950SBarry Smith   ierr  = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
4218965ea79SLois Curfman McInnes   count = 0;
4228965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4238965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4248965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4258965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4268965ea79SLois Curfman McInnes     }
4278965ea79SLois Curfman McInnes   }
428606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
42974ed9c26SBarry Smith   ierr = PetscFree2(lens,source);CHKERRQ(ierr);
430606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
431606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
4328965ea79SLois Curfman McInnes 
43397b48c8fSBarry Smith   /* fix right hand side if needed */
43497b48c8fSBarry Smith   if (x && b) {
43597b48c8fSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
43697b48c8fSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
43797b48c8fSBarry Smith     for (i=0; i<slen; i++) {
43897b48c8fSBarry Smith       bb[lrows[i]] = diag*xx[lrows[i]];
43997b48c8fSBarry Smith     }
44097b48c8fSBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
44197b48c8fSBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
44297b48c8fSBarry Smith   }
44397b48c8fSBarry Smith 
4448965ea79SLois Curfman McInnes   /* actually zap the local rows */
445b9679d65SBarry Smith   ierr = MatZeroRows(l->A,slen,lrows,0.0,0,0);CHKERRQ(ierr);
446e2eb51b1SBarry Smith   if (diag != 0.0) {
447b9679d65SBarry Smith     Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data;
448b9679d65SBarry Smith     PetscInt     m   = ll->lda, i;
449b9679d65SBarry Smith 
450b9679d65SBarry Smith     for (i=0; i<slen; i++) {
451b9679d65SBarry Smith       ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag;
452b9679d65SBarry Smith     }
453b9679d65SBarry Smith   }
454606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4558965ea79SLois Curfman McInnes 
4568965ea79SLois Curfman McInnes   /* wait on sends */
4578965ea79SLois Curfman McInnes   if (nsends) {
458b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
459ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
460606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4618965ea79SLois Curfman McInnes   }
462606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
463606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4643a40ed3dSBarry Smith   PetscFunctionReturn(0);
4658965ea79SLois Curfman McInnes }
4668965ea79SLois Curfman McInnes 
4674a2ae208SSatish Balay #undef __FUNCT__
4684a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
469dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4708965ea79SLois Curfman McInnes {
47139ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
472dfbe8321SBarry Smith   PetscErrorCode ierr;
473c456f294SBarry Smith 
4743a40ed3dSBarry Smith   PetscFunctionBegin;
475ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
476ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
47744cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4783a40ed3dSBarry Smith   PetscFunctionReturn(0);
4798965ea79SLois Curfman McInnes }
4808965ea79SLois Curfman McInnes 
4814a2ae208SSatish Balay #undef __FUNCT__
4824a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
483dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4848965ea79SLois Curfman McInnes {
48539ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
486dfbe8321SBarry Smith   PetscErrorCode ierr;
487c456f294SBarry Smith 
4883a40ed3dSBarry Smith   PetscFunctionBegin;
489ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
490ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
49144cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4923a40ed3dSBarry Smith   PetscFunctionReturn(0);
4938965ea79SLois Curfman McInnes }
4948965ea79SLois Curfman McInnes 
4954a2ae208SSatish Balay #undef __FUNCT__
4964a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
497dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
498096963f5SLois Curfman McInnes {
499096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
500dfbe8321SBarry Smith   PetscErrorCode ierr;
50187828ca2SBarry Smith   PetscScalar    zero = 0.0;
502096963f5SLois Curfman McInnes 
5033a40ed3dSBarry Smith   PetscFunctionBegin;
5042dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
5057c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
506ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
507ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5083a40ed3dSBarry Smith   PetscFunctionReturn(0);
509096963f5SLois Curfman McInnes }
510096963f5SLois Curfman McInnes 
5114a2ae208SSatish Balay #undef __FUNCT__
5124a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
513dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
514096963f5SLois Curfman McInnes {
515096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
516dfbe8321SBarry Smith   PetscErrorCode ierr;
517096963f5SLois Curfman McInnes 
5183a40ed3dSBarry Smith   PetscFunctionBegin;
5193501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
5207c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
521ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
522ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5233a40ed3dSBarry Smith   PetscFunctionReturn(0);
524096963f5SLois Curfman McInnes }
525096963f5SLois Curfman McInnes 
5264a2ae208SSatish Balay #undef __FUNCT__
5274a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
528dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5298965ea79SLois Curfman McInnes {
53039ddd567SLois Curfman McInnes   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
531096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
532dfbe8321SBarry Smith   PetscErrorCode ierr;
533d0f46423SBarry Smith   PetscInt       len,i,n,m = A->rmap->n,radd;
53487828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
535ed3cc1f0SBarry Smith 
5363a40ed3dSBarry Smith   PetscFunctionBegin;
5372dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5381ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
539096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
540e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
541d0f46423SBarry Smith   len  = PetscMin(a->A->rmap->n,a->A->cmap->n);
542d0f46423SBarry Smith   radd = A->rmap->rstart*m;
54344cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
544096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
545096963f5SLois Curfman McInnes   }
5461ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5473a40ed3dSBarry Smith   PetscFunctionReturn(0);
5488965ea79SLois Curfman McInnes }
5498965ea79SLois Curfman McInnes 
5504a2ae208SSatish Balay #undef __FUNCT__
5514a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
552dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5538965ea79SLois Curfman McInnes {
5543501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
555dfbe8321SBarry Smith   PetscErrorCode ierr;
556ed3cc1f0SBarry Smith 
5573a40ed3dSBarry Smith   PetscFunctionBegin;
558aa482453SBarry Smith #if defined(PETSC_USE_LOG)
559d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N);
5608965ea79SLois Curfman McInnes #endif
5618798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
5626bf464f9SBarry Smith   ierr = MatDestroy(&mdn->A);CHKERRQ(ierr);
5636bf464f9SBarry Smith   ierr = VecDestroy(&mdn->lvec);CHKERRQ(ierr);
5646bf464f9SBarry Smith   ierr = VecScatterDestroy(&mdn->Mvctx);CHKERRQ(ierr);
56501b82886SBarry Smith 
566bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
567dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
5680298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",NULL);CHKERRQ(ierr);
5690298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",NULL);CHKERRQ(ierr);
5700298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","",NULL);CHKERRQ(ierr);
5710298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","",NULL);CHKERRQ(ierr);
5720298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","",NULL);CHKERRQ(ierr);
5733a40ed3dSBarry Smith   PetscFunctionReturn(0);
5748965ea79SLois Curfman McInnes }
57539ddd567SLois Curfman McInnes 
5764a2ae208SSatish Balay #undef __FUNCT__
5774a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5786849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5798965ea79SLois Curfman McInnes {
58039ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
581dfbe8321SBarry Smith   PetscErrorCode    ierr;
582aa05aa95SBarry Smith   PetscViewerFormat format;
583aa05aa95SBarry Smith   int               fd;
584d0f46423SBarry Smith   PetscInt          header[4],mmax,N = mat->cmap->N,i,j,m,k;
585aa05aa95SBarry Smith   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
586578230a0SSatish Balay   PetscScalar       *work,*v,*vv;
587aa05aa95SBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;
5887056b6fcSBarry Smith 
5893a40ed3dSBarry Smith   PetscFunctionBegin;
59039ddd567SLois Curfman McInnes   if (mdn->size == 1) {
59139ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
592aa05aa95SBarry Smith   } else {
593aa05aa95SBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
594*ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr);
595*ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
596aa05aa95SBarry Smith 
597aa05aa95SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
598f4403165SShri Abhyankar     if (format == PETSC_VIEWER_NATIVE) {
599aa05aa95SBarry Smith 
600aa05aa95SBarry Smith       if (!rank) {
601aa05aa95SBarry Smith         /* store the matrix as a dense matrix */
6020700a824SBarry Smith         header[0] = MAT_FILE_CLASSID;
603d0f46423SBarry Smith         header[1] = mat->rmap->N;
604aa05aa95SBarry Smith         header[2] = N;
605aa05aa95SBarry Smith         header[3] = MATRIX_BINARY_FORMAT_DENSE;
606aa05aa95SBarry Smith         ierr      = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
607aa05aa95SBarry Smith 
608aa05aa95SBarry Smith         /* get largest work array needed for transposing array */
609d0f46423SBarry Smith         mmax = mat->rmap->n;
610aa05aa95SBarry Smith         for (i=1; i<size; i++) {
611d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
6128965ea79SLois Curfman McInnes         }
613aa05aa95SBarry Smith         ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&work);CHKERRQ(ierr);
614aa05aa95SBarry Smith 
615aa05aa95SBarry Smith         /* write out local array, by rows */
616d0f46423SBarry Smith         m = mat->rmap->n;
617aa05aa95SBarry Smith         v = a->v;
618aa05aa95SBarry Smith         for (j=0; j<N; j++) {
619aa05aa95SBarry Smith           for (i=0; i<m; i++) {
620578230a0SSatish Balay             work[j + i*N] = *v++;
621aa05aa95SBarry Smith           }
622aa05aa95SBarry Smith         }
623aa05aa95SBarry Smith         ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
624aa05aa95SBarry Smith         /* get largest work array to receive messages from other processes, excludes process zero */
625aa05aa95SBarry Smith         mmax = 0;
626aa05aa95SBarry Smith         for (i=1; i<size; i++) {
627d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
628aa05aa95SBarry Smith         }
629578230a0SSatish Balay         ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&vv);CHKERRQ(ierr);
630aa05aa95SBarry Smith         for (k = 1; k < size; k++) {
631f8009846SMatthew Knepley           v    = vv;
632d0f46423SBarry Smith           m    = mat->rmap->range[k+1] - mat->rmap->range[k];
633*ce94432eSBarry Smith           ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
634aa05aa95SBarry Smith 
635aa05aa95SBarry Smith           for (j = 0; j < N; j++) {
636aa05aa95SBarry Smith             for (i = 0; i < m; i++) {
637578230a0SSatish Balay               work[j + i*N] = *v++;
638aa05aa95SBarry Smith             }
639aa05aa95SBarry Smith           }
640aa05aa95SBarry Smith           ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
641aa05aa95SBarry Smith         }
642aa05aa95SBarry Smith         ierr = PetscFree(work);CHKERRQ(ierr);
643578230a0SSatish Balay         ierr = PetscFree(vv);CHKERRQ(ierr);
644aa05aa95SBarry Smith       } else {
645*ce94432eSBarry Smith         ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
646aa05aa95SBarry Smith       }
647eb3f98d2SBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE)");
648aa05aa95SBarry Smith   }
6493a40ed3dSBarry Smith   PetscFunctionReturn(0);
6508965ea79SLois Curfman McInnes }
6518965ea79SLois Curfman McInnes 
6524a2ae208SSatish Balay #undef __FUNCT__
6534a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
6546849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6558965ea79SLois Curfman McInnes {
65639ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
657dfbe8321SBarry Smith   PetscErrorCode    ierr;
65832dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
65919fd82e9SBarry Smith   PetscViewerType   vtype;
660ace3abfcSBarry Smith   PetscBool         iascii,isdraw;
661b0a32e0cSBarry Smith   PetscViewer       sviewer;
662f3ef73ceSBarry Smith   PetscViewerFormat format;
6638965ea79SLois Curfman McInnes 
6643a40ed3dSBarry Smith   PetscFunctionBegin;
665251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
666251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
66732077d6dSBarry Smith   if (iascii) {
668b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
669b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
670456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6714e220ebcSLois Curfman McInnes       MatInfo info;
672888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
6737b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
6747b23a99aSBarry 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);
675b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6767b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
6775d00a290SHong Zhang       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
6783a40ed3dSBarry Smith       PetscFunctionReturn(0);
679fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6803a40ed3dSBarry Smith       PetscFunctionReturn(0);
6818965ea79SLois Curfman McInnes     }
682f1af5d2fSBarry Smith   } else if (isdraw) {
683b0a32e0cSBarry Smith     PetscDraw draw;
684ace3abfcSBarry Smith     PetscBool isnull;
685f1af5d2fSBarry Smith 
686b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
687b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
688f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
689f1af5d2fSBarry Smith   }
69077ed5343SBarry Smith 
6918965ea79SLois Curfman McInnes   if (size == 1) {
69239ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6933a40ed3dSBarry Smith   } else {
6948965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
6958965ea79SLois Curfman McInnes     Mat         A;
696d0f46423SBarry Smith     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
697ba8c8a56SBarry Smith     PetscInt    *cols;
698ba8c8a56SBarry Smith     PetscScalar *vals;
6998965ea79SLois Curfman McInnes 
700*ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr);
7018965ea79SLois Curfman McInnes     if (!rank) {
702f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
7033a40ed3dSBarry Smith     } else {
704f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
7058965ea79SLois Curfman McInnes     }
7067adad957SLisandro Dalcin     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
707878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
7080298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(A,NULL);CHKERRQ(ierr);
70952e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
7108965ea79SLois Curfman McInnes 
71139ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
71239ddd567SLois Curfman McInnes        but it's quick for now */
71351022da4SBarry Smith     A->insertmode = INSERT_VALUES;
7142205254eSKarl Rupp 
7152205254eSKarl Rupp     row = mat->rmap->rstart;
7162205254eSKarl Rupp     m   = mdn->A->rmap->n;
7178965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
718ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
719ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
720ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
72139ddd567SLois Curfman McInnes       row++;
7228965ea79SLois Curfman McInnes     }
7238965ea79SLois Curfman McInnes 
7246d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7256d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
726b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
727b9b97703SBarry Smith     if (!rank) {
7287566de4bSShri Abhyankar       ierr = PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr);
7297566de4bSShri Abhyankar       /* Set the type name to MATMPIDense so that the correct type can be printed out by PetscObjectPrintClassNamePrefixType() in MatView_SeqDense_ASCII()*/
7307566de4bSShri Abhyankar       PetscStrcpy(((PetscObject)((Mat_MPIDense*)(A->data))->A)->type_name,MATMPIDENSE);
7316831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
7328965ea79SLois Curfman McInnes     }
733b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
734b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7356bf464f9SBarry Smith     ierr = MatDestroy(&A);CHKERRQ(ierr);
7368965ea79SLois Curfman McInnes   }
7373a40ed3dSBarry Smith   PetscFunctionReturn(0);
7388965ea79SLois Curfman McInnes }
7398965ea79SLois Curfman McInnes 
7404a2ae208SSatish Balay #undef __FUNCT__
7414a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
742dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
7438965ea79SLois Curfman McInnes {
744dfbe8321SBarry Smith   PetscErrorCode ierr;
745ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw,issocket;
7468965ea79SLois Curfman McInnes 
747433994e6SBarry Smith   PetscFunctionBegin;
748251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
749251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
750251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr);
751251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
7520f5bd95cSBarry Smith 
75332077d6dSBarry Smith   if (iascii || issocket || isdraw) {
754f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7550f5bd95cSBarry Smith   } else if (isbinary) {
7563a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
75711aeaf0aSBarry Smith   }
7583a40ed3dSBarry Smith   PetscFunctionReturn(0);
7598965ea79SLois Curfman McInnes }
7608965ea79SLois Curfman McInnes 
7614a2ae208SSatish Balay #undef __FUNCT__
7624a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
763dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7648965ea79SLois Curfman McInnes {
7653501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7663501a2bdSLois Curfman McInnes   Mat            mdn  = mat->A;
767dfbe8321SBarry Smith   PetscErrorCode ierr;
768329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7698965ea79SLois Curfman McInnes 
7703a40ed3dSBarry Smith   PetscFunctionBegin;
7714e220ebcSLois Curfman McInnes   info->block_size = 1.0;
7722205254eSKarl Rupp 
7734e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7742205254eSKarl Rupp 
7754e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7764e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7778965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7784e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7794e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7804e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7814e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7824e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7838965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
784*ce94432eSBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
7852205254eSKarl Rupp 
7864e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7874e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7884e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7894e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7904e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7918965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
792*ce94432eSBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
7932205254eSKarl Rupp 
7944e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7954e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7964e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7974e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7984e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7998965ea79SLois Curfman McInnes   }
8004e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
8014e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
8024e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
8033a40ed3dSBarry Smith   PetscFunctionReturn(0);
8048965ea79SLois Curfman McInnes }
8058965ea79SLois Curfman McInnes 
8064a2ae208SSatish Balay #undef __FUNCT__
8074a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
808ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg)
8098965ea79SLois Curfman McInnes {
81039ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
811dfbe8321SBarry Smith   PetscErrorCode ierr;
8128965ea79SLois Curfman McInnes 
8133a40ed3dSBarry Smith   PetscFunctionBegin;
81412c028f9SKris Buschelman   switch (op) {
815512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
81612c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
81712c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
8184e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
81912c028f9SKris Buschelman     break;
82012c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
8214e0d8c25SBarry Smith     a->roworiented = flg;
8222205254eSKarl Rupp 
8234e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82412c028f9SKris Buschelman     break;
8254e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
82613fa8e87SLisandro Dalcin   case MAT_KEEP_NONZERO_PATTERN:
82712c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
828290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
82912c028f9SKris Buschelman     break;
83012c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
8314e0d8c25SBarry Smith     a->donotstash = flg;
83212c028f9SKris Buschelman     break;
83377e54ba9SKris Buschelman   case MAT_SYMMETRIC:
83477e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8359a4540c5SBarry Smith   case MAT_HERMITIAN:
8369a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
837600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
838290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
83977e54ba9SKris Buschelman     break;
84012c028f9SKris Buschelman   default:
841e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
8423a40ed3dSBarry Smith   }
8433a40ed3dSBarry Smith   PetscFunctionReturn(0);
8448965ea79SLois Curfman McInnes }
8458965ea79SLois Curfman McInnes 
8468965ea79SLois Curfman McInnes 
8474a2ae208SSatish Balay #undef __FUNCT__
8484a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
849dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8505b2fa520SLois Curfman McInnes {
8515b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8525b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
85387828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
854dfbe8321SBarry Smith   PetscErrorCode ierr;
855d0f46423SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n;
8565b2fa520SLois Curfman McInnes 
8575b2fa520SLois Curfman McInnes   PetscFunctionBegin;
85872d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8595b2fa520SLois Curfman McInnes   if (ll) {
86072d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
861e32f2f54SBarry Smith     if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8621ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8635b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8645b2fa520SLois Curfman McInnes       x = l[i];
8655b2fa520SLois Curfman McInnes       v = mat->v + i;
8665b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8675b2fa520SLois Curfman McInnes     }
8681ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
869efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8705b2fa520SLois Curfman McInnes   }
8715b2fa520SLois Curfman McInnes   if (rr) {
872175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
873e32f2f54SBarry Smith     if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
874ca9f406cSSatish Balay     ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
875ca9f406cSSatish Balay     ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8761ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8775b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8785b2fa520SLois Curfman McInnes       x = r[i];
8795b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8802205254eSKarl Rupp       for (j=0; j<m; j++) (*v++) *= x;
8815b2fa520SLois Curfman McInnes     }
8821ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
883efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8845b2fa520SLois Curfman McInnes   }
8855b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8865b2fa520SLois Curfman McInnes }
8875b2fa520SLois Curfman McInnes 
8884a2ae208SSatish Balay #undef __FUNCT__
8894a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
890dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
891096963f5SLois Curfman McInnes {
8923501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8933501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
894dfbe8321SBarry Smith   PetscErrorCode ierr;
89513f74950SBarry Smith   PetscInt       i,j;
896329f5518SBarry Smith   PetscReal      sum = 0.0;
89787828ca2SBarry Smith   PetscScalar    *v  = mat->v;
8983501a2bdSLois Curfman McInnes 
8993a40ed3dSBarry Smith   PetscFunctionBegin;
9003501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
901064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
9023501a2bdSLois Curfman McInnes   } else {
9033501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
904d0f46423SBarry Smith       for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) {
905329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
9063501a2bdSLois Curfman McInnes       }
907*ce94432eSBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
9088f1a2a5eSBarry Smith       *nrm = PetscSqrtReal(*nrm);
909dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr);
9103a40ed3dSBarry Smith     } else if (type == NORM_1) {
911329f5518SBarry Smith       PetscReal *tmp,*tmp2;
91274ed9c26SBarry Smith       ierr = PetscMalloc2(A->cmap->N,PetscReal,&tmp,A->cmap->N,PetscReal,&tmp2);CHKERRQ(ierr);
91374ed9c26SBarry Smith       ierr = PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
91474ed9c26SBarry Smith       ierr = PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
915064f8208SBarry Smith       *nrm = 0.0;
9163501a2bdSLois Curfman McInnes       v    = mat->v;
917d0f46423SBarry Smith       for (j=0; j<mdn->A->cmap->n; j++) {
918d0f46423SBarry Smith         for (i=0; i<mdn->A->rmap->n; i++) {
91967e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
9203501a2bdSLois Curfman McInnes         }
9213501a2bdSLois Curfman McInnes       }
922*ce94432eSBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
923d0f46423SBarry Smith       for (j=0; j<A->cmap->N; j++) {
924064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
9253501a2bdSLois Curfman McInnes       }
92674ed9c26SBarry Smith       ierr = PetscFree2(tmp,tmp);CHKERRQ(ierr);
927d0f46423SBarry Smith       ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
9283a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
929329f5518SBarry Smith       PetscReal ntemp;
9303501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
931*ce94432eSBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
932*ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm");
9333501a2bdSLois Curfman McInnes   }
9343a40ed3dSBarry Smith   PetscFunctionReturn(0);
9353501a2bdSLois Curfman McInnes }
9363501a2bdSLois Curfman McInnes 
9374a2ae208SSatish Balay #undef __FUNCT__
9384a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
939fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
9403501a2bdSLois Curfman McInnes {
9413501a2bdSLois Curfman McInnes   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
9423501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
9433501a2bdSLois Curfman McInnes   Mat            B;
944d0f46423SBarry Smith   PetscInt       M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart;
9456849ba73SBarry Smith   PetscErrorCode ierr;
94613f74950SBarry Smith   PetscInt       j,i;
94787828ca2SBarry Smith   PetscScalar    *v;
9483501a2bdSLois Curfman McInnes 
9493a40ed3dSBarry Smith   PetscFunctionBegin;
950*ce94432eSBarry Smith   if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place");
951fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
952*ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
953d0f46423SBarry Smith     ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr);
9547adad957SLisandro Dalcin     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
9550298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr);
956fc4dec0aSBarry Smith   } else {
957fc4dec0aSBarry Smith     B = *matout;
958fc4dec0aSBarry Smith   }
9593501a2bdSLois Curfman McInnes 
960d0f46423SBarry Smith   m    = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v;
9611acff37aSSatish Balay   ierr = PetscMalloc(m*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
9623501a2bdSLois Curfman McInnes   for (i=0; i<m; i++) rwork[i] = rstart + i;
9631acff37aSSatish Balay   for (j=0; j<n; j++) {
9643501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9653501a2bdSLois Curfman McInnes     v   += m;
9663501a2bdSLois Curfman McInnes   }
967606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9686d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9696d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
970815cbec1SBarry Smith   if (reuse == MAT_INITIAL_MATRIX || *matout != A) {
9713501a2bdSLois Curfman McInnes     *matout = B;
9723501a2bdSLois Curfman McInnes   } else {
973eb6b5d47SBarry Smith     ierr = MatHeaderMerge(A,B);CHKERRQ(ierr);
9743501a2bdSLois Curfman McInnes   }
9753a40ed3dSBarry Smith   PetscFunctionReturn(0);
976096963f5SLois Curfman McInnes }
977096963f5SLois Curfman McInnes 
97844cd7ae7SLois Curfman McInnes 
9796849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*);
980d84338a6SBarry Smith extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar);
9818965ea79SLois Curfman McInnes 
9824a2ae208SSatish Balay #undef __FUNCT__
9834994cf47SJed Brown #define __FUNCT__ "MatSetUp_MPIDense"
9844994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A)
985273d9f13SBarry Smith {
986dfbe8321SBarry Smith   PetscErrorCode ierr;
987273d9f13SBarry Smith 
988273d9f13SBarry Smith   PetscFunctionBegin;
989273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
990273d9f13SBarry Smith   PetscFunctionReturn(0);
991273d9f13SBarry Smith }
992273d9f13SBarry Smith 
99330716080SHong Zhang #undef __FUNCT__
994488007eeSBarry Smith #define __FUNCT__ "MatAXPY_MPIDense"
995488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
996488007eeSBarry Smith {
997488007eeSBarry Smith   PetscErrorCode ierr;
998488007eeSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;
999488007eeSBarry Smith 
1000488007eeSBarry Smith   PetscFunctionBegin;
1001488007eeSBarry Smith   ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr);
1002488007eeSBarry Smith   PetscFunctionReturn(0);
1003488007eeSBarry Smith }
1004488007eeSBarry Smith 
1005ba337c44SJed Brown #undef __FUNCT__
1006ba337c44SJed Brown #define __FUNCT__ "MatConjugate_MPIDense"
10077087cfbeSBarry Smith PetscErrorCode  MatConjugate_MPIDense(Mat mat)
1008ba337c44SJed Brown {
1009ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)mat->data;
1010ba337c44SJed Brown   PetscErrorCode ierr;
1011ba337c44SJed Brown 
1012ba337c44SJed Brown   PetscFunctionBegin;
1013ba337c44SJed Brown   ierr = MatConjugate(a->A);CHKERRQ(ierr);
1014ba337c44SJed Brown   PetscFunctionReturn(0);
1015ba337c44SJed Brown }
1016ba337c44SJed Brown 
1017ba337c44SJed Brown #undef __FUNCT__
1018ba337c44SJed Brown #define __FUNCT__ "MatRealPart_MPIDense"
1019ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A)
1020ba337c44SJed Brown {
1021ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1022ba337c44SJed Brown   PetscErrorCode ierr;
1023ba337c44SJed Brown 
1024ba337c44SJed Brown   PetscFunctionBegin;
1025ba337c44SJed Brown   ierr = MatRealPart(a->A);CHKERRQ(ierr);
1026ba337c44SJed Brown   PetscFunctionReturn(0);
1027ba337c44SJed Brown }
1028ba337c44SJed Brown 
1029ba337c44SJed Brown #undef __FUNCT__
1030ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_MPIDense"
1031ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1032ba337c44SJed Brown {
1033ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1034ba337c44SJed Brown   PetscErrorCode ierr;
1035ba337c44SJed Brown 
1036ba337c44SJed Brown   PetscFunctionBegin;
1037ba337c44SJed Brown   ierr = MatImaginaryPart(a->A);CHKERRQ(ierr);
1038ba337c44SJed Brown   PetscFunctionReturn(0);
1039ba337c44SJed Brown }
1040ba337c44SJed Brown 
10410716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);
10420716a85fSBarry Smith #undef __FUNCT__
10430716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense"
10440716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
10450716a85fSBarry Smith {
10460716a85fSBarry Smith   PetscErrorCode ierr;
10470716a85fSBarry Smith   PetscInt       i,n;
10480716a85fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
10490716a85fSBarry Smith   PetscReal      *work;
10500716a85fSBarry Smith 
10510716a85fSBarry Smith   PetscFunctionBegin;
10520298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr);
10530716a85fSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr);
10540716a85fSBarry Smith   ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr);
10550716a85fSBarry Smith   if (type == NORM_2) {
10560716a85fSBarry Smith     for (i=0; i<n; i++) work[i] *= work[i];
10570716a85fSBarry Smith   }
10580716a85fSBarry Smith   if (type == NORM_INFINITY) {
10590716a85fSBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr);
10600716a85fSBarry Smith   } else {
10610716a85fSBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr);
10620716a85fSBarry Smith   }
10630716a85fSBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
10640716a85fSBarry Smith   if (type == NORM_2) {
10658f1a2a5eSBarry Smith     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
10660716a85fSBarry Smith   }
10670716a85fSBarry Smith   PetscFunctionReturn(0);
10680716a85fSBarry Smith }
10690716a85fSBarry Smith 
107073a71a0fSBarry Smith #undef __FUNCT__
107173a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_MPIDense"
107273a71a0fSBarry Smith static PetscErrorCode  MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
107373a71a0fSBarry Smith {
107473a71a0fSBarry Smith   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;
107573a71a0fSBarry Smith   PetscErrorCode ierr;
107673a71a0fSBarry Smith   PetscScalar    *a;
107773a71a0fSBarry Smith   PetscInt       m,n,i;
107873a71a0fSBarry Smith 
107973a71a0fSBarry Smith   PetscFunctionBegin;
108073a71a0fSBarry Smith   ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr);
10818c778c55SBarry Smith   ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr);
108273a71a0fSBarry Smith   for (i=0; i<m*n; i++) {
108373a71a0fSBarry Smith     ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr);
108473a71a0fSBarry Smith   }
10858c778c55SBarry Smith   ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr);
108673a71a0fSBarry Smith   PetscFunctionReturn(0);
108773a71a0fSBarry Smith }
108873a71a0fSBarry Smith 
10898965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
109009dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense,
109109dc0095SBarry Smith                                         MatGetRow_MPIDense,
109209dc0095SBarry Smith                                         MatRestoreRow_MPIDense,
109309dc0095SBarry Smith                                         MatMult_MPIDense,
109497304618SKris Buschelman                                 /*  4*/ MatMultAdd_MPIDense,
10957c922b88SBarry Smith                                         MatMultTranspose_MPIDense,
10967c922b88SBarry Smith                                         MatMultTransposeAdd_MPIDense,
10978965ea79SLois Curfman McInnes                                         0,
109809dc0095SBarry Smith                                         0,
109909dc0095SBarry Smith                                         0,
110097304618SKris Buschelman                                 /* 10*/ 0,
110109dc0095SBarry Smith                                         0,
110209dc0095SBarry Smith                                         0,
110309dc0095SBarry Smith                                         0,
110409dc0095SBarry Smith                                         MatTranspose_MPIDense,
110597304618SKris Buschelman                                 /* 15*/ MatGetInfo_MPIDense,
11066e4ee0c6SHong Zhang                                         MatEqual_MPIDense,
110709dc0095SBarry Smith                                         MatGetDiagonal_MPIDense,
11085b2fa520SLois Curfman McInnes                                         MatDiagonalScale_MPIDense,
110909dc0095SBarry Smith                                         MatNorm_MPIDense,
111097304618SKris Buschelman                                 /* 20*/ MatAssemblyBegin_MPIDense,
111109dc0095SBarry Smith                                         MatAssemblyEnd_MPIDense,
111209dc0095SBarry Smith                                         MatSetOption_MPIDense,
111309dc0095SBarry Smith                                         MatZeroEntries_MPIDense,
1114d519adbfSMatthew Knepley                                 /* 24*/ MatZeroRows_MPIDense,
1115919b68f7SBarry Smith                                         0,
111601b82886SBarry Smith                                         0,
111701b82886SBarry Smith                                         0,
111801b82886SBarry Smith                                         0,
11194994cf47SJed Brown                                 /* 29*/ MatSetUp_MPIDense,
1120273d9f13SBarry Smith                                         0,
112109dc0095SBarry Smith                                         0,
11228c778c55SBarry Smith                                         0,
11238c778c55SBarry Smith                                         0,
1124d519adbfSMatthew Knepley                                 /* 34*/ MatDuplicate_MPIDense,
112509dc0095SBarry Smith                                         0,
112609dc0095SBarry Smith                                         0,
112709dc0095SBarry Smith                                         0,
112809dc0095SBarry Smith                                         0,
1129d519adbfSMatthew Knepley                                 /* 39*/ MatAXPY_MPIDense,
11302ce60cd0SSatish Balay                                         MatGetSubMatrices_MPIDense,
113109dc0095SBarry Smith                                         0,
113209dc0095SBarry Smith                                         MatGetValues_MPIDense,
113309dc0095SBarry Smith                                         0,
1134d519adbfSMatthew Knepley                                 /* 44*/ 0,
113509dc0095SBarry Smith                                         MatScale_MPIDense,
113609dc0095SBarry Smith                                         0,
113709dc0095SBarry Smith                                         0,
113809dc0095SBarry Smith                                         0,
113973a71a0fSBarry Smith                                 /* 49*/ MatSetRandom_MPIDense,
114009dc0095SBarry Smith                                         0,
114109dc0095SBarry Smith                                         0,
114209dc0095SBarry Smith                                         0,
114309dc0095SBarry Smith                                         0,
1144d519adbfSMatthew Knepley                                 /* 54*/ 0,
114509dc0095SBarry Smith                                         0,
114609dc0095SBarry Smith                                         0,
114709dc0095SBarry Smith                                         0,
114809dc0095SBarry Smith                                         0,
1149d519adbfSMatthew Knepley                                 /* 59*/ MatGetSubMatrix_MPIDense,
1150b9b97703SBarry Smith                                         MatDestroy_MPIDense,
1151b9b97703SBarry Smith                                         MatView_MPIDense,
1152357abbc8SBarry Smith                                         0,
115397304618SKris Buschelman                                         0,
1154d519adbfSMatthew Knepley                                 /* 64*/ 0,
115597304618SKris Buschelman                                         0,
115697304618SKris Buschelman                                         0,
115797304618SKris Buschelman                                         0,
115897304618SKris Buschelman                                         0,
1159d519adbfSMatthew Knepley                                 /* 69*/ 0,
116097304618SKris Buschelman                                         0,
116197304618SKris Buschelman                                         0,
116297304618SKris Buschelman                                         0,
116397304618SKris Buschelman                                         0,
1164d519adbfSMatthew Knepley                                 /* 74*/ 0,
116597304618SKris Buschelman                                         0,
116697304618SKris Buschelman                                         0,
116797304618SKris Buschelman                                         0,
116897304618SKris Buschelman                                         0,
1169d519adbfSMatthew Knepley                                 /* 79*/ 0,
117097304618SKris Buschelman                                         0,
117197304618SKris Buschelman                                         0,
117297304618SKris Buschelman                                         0,
11735bba2384SShri Abhyankar                                 /* 83*/ MatLoad_MPIDense,
1174865e5f61SKris Buschelman                                         0,
1175865e5f61SKris Buschelman                                         0,
1176865e5f61SKris Buschelman                                         0,
1177865e5f61SKris Buschelman                                         0,
1178865e5f61SKris Buschelman                                         0,
1179d519adbfSMatthew Knepley                                 /* 89*/
1180865e5f61SKris Buschelman                                         0,
1181865e5f61SKris Buschelman                                         0,
1182865e5f61SKris Buschelman                                         0,
11832fbe02b9SBarry Smith                                         0,
1184ba337c44SJed Brown                                         0,
1185d519adbfSMatthew Knepley                                 /* 94*/ 0,
1186865e5f61SKris Buschelman                                         0,
1187865e5f61SKris Buschelman                                         0,
1188ba337c44SJed Brown                                         0,
1189ba337c44SJed Brown                                         0,
1190ba337c44SJed Brown                                 /* 99*/ 0,
1191ba337c44SJed Brown                                         0,
1192ba337c44SJed Brown                                         0,
1193ba337c44SJed Brown                                         MatConjugate_MPIDense,
1194ba337c44SJed Brown                                         0,
1195ba337c44SJed Brown                                 /*104*/ 0,
1196ba337c44SJed Brown                                         MatRealPart_MPIDense,
1197ba337c44SJed Brown                                         MatImaginaryPart_MPIDense,
119886d161a7SShri Abhyankar                                         0,
119986d161a7SShri Abhyankar                                         0,
120086d161a7SShri Abhyankar                                 /*109*/ 0,
120186d161a7SShri Abhyankar                                         0,
120286d161a7SShri Abhyankar                                         0,
120386d161a7SShri Abhyankar                                         0,
120486d161a7SShri Abhyankar                                         0,
120586d161a7SShri Abhyankar                                 /*114*/ 0,
120686d161a7SShri Abhyankar                                         0,
120786d161a7SShri Abhyankar                                         0,
120886d161a7SShri Abhyankar                                         0,
120986d161a7SShri Abhyankar                                         0,
121086d161a7SShri Abhyankar                                 /*119*/ 0,
121186d161a7SShri Abhyankar                                         0,
121286d161a7SShri Abhyankar                                         0,
12130716a85fSBarry Smith                                         0,
12140716a85fSBarry Smith                                         0,
12150716a85fSBarry Smith                                 /*124*/ 0,
12160716a85fSBarry Smith                                         MatGetColumnNorms_MPIDense
1217ba337c44SJed Brown };
12188965ea79SLois Curfman McInnes 
1219273d9f13SBarry Smith EXTERN_C_BEGIN
12204a2ae208SSatish Balay #undef __FUNCT__
1221a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
12227087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1223a23d5eceSKris Buschelman {
1224a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1225dfbe8321SBarry Smith   PetscErrorCode ierr;
1226a23d5eceSKris Buschelman 
1227a23d5eceSKris Buschelman   PetscFunctionBegin;
1228a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1229a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1230a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1231a23d5eceSKris Buschelman 
1232a23d5eceSKris Buschelman   a       = (Mat_MPIDense*)mat->data;
123334ef9618SShri Abhyankar   ierr    = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
123434ef9618SShri Abhyankar   ierr    = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
123534ef9618SShri Abhyankar   a->nvec = mat->cmap->n;
123634ef9618SShri Abhyankar 
1237f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1238d0f46423SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr);
12395c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
12405c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
124152e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1242a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1243a23d5eceSKris Buschelman }
1244a23d5eceSKris Buschelman EXTERN_C_END
1245a23d5eceSKris Buschelman 
1246a23d5eceSKris Buschelman EXTERN_C_BEGIN
1247a23d5eceSKris Buschelman #undef __FUNCT__
12484a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
12497087cfbeSBarry Smith PetscErrorCode  MatCreate_MPIDense(Mat mat)
1250273d9f13SBarry Smith {
1251273d9f13SBarry Smith   Mat_MPIDense   *a;
1252dfbe8321SBarry Smith   PetscErrorCode ierr;
1253273d9f13SBarry Smith 
1254273d9f13SBarry Smith   PetscFunctionBegin;
125538f2d2fdSLisandro Dalcin   ierr      = PetscNewLog(mat,Mat_MPIDense,&a);CHKERRQ(ierr);
1256b0a32e0cSBarry Smith   mat->data = (void*)a;
1257273d9f13SBarry Smith   ierr      = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1258273d9f13SBarry Smith 
1259273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1260*ce94432eSBarry Smith   ierr            = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);CHKERRQ(ierr);
1261*ce94432eSBarry Smith   ierr            = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);CHKERRQ(ierr);
1262273d9f13SBarry Smith 
1263273d9f13SBarry Smith   /* build cache for off array entries formed */
1264273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
12652205254eSKarl Rupp 
1266*ce94432eSBarry Smith   ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);CHKERRQ(ierr);
1267273d9f13SBarry Smith 
1268273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1269273d9f13SBarry Smith   a->lvec        = 0;
1270273d9f13SBarry Smith   a->Mvctx       = 0;
1271273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1272273d9f13SBarry Smith 
12738c778c55SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseGetArray_C","MatDenseGetArray_MPIDense",MatDenseGetArray_MPIDense);CHKERRQ(ierr);
12748c778c55SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseRestoreArray_C","MatDenseRestoreArray_MPIDense",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr);
12758c778c55SBarry Smith 
1276273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1277273d9f13SBarry Smith                                            "MatGetDiagonalBlock_MPIDense",
1278273d9f13SBarry Smith                                            MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1279a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1280a23d5eceSKris Buschelman                                            "MatMPIDenseSetPreallocation_MPIDense",
1281a23d5eceSKris Buschelman                                            MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
12824ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",
12834ae313f4SHong Zhang                                            "MatMatMult_MPIAIJ_MPIDense",
12844ae313f4SHong Zhang                                            MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
12854ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",
12864ae313f4SHong Zhang                                            "MatMatMultSymbolic_MPIAIJ_MPIDense",
12874ae313f4SHong Zhang                                            MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
12884ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",
12894ae313f4SHong Zhang                                            "MatMatMultNumeric_MPIAIJ_MPIDense",
12904ae313f4SHong Zhang                                            MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
129138aed534SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr);
1292273d9f13SBarry Smith   PetscFunctionReturn(0);
1293273d9f13SBarry Smith }
1294273d9f13SBarry Smith EXTERN_C_END
1295273d9f13SBarry Smith 
1296209238afSKris Buschelman /*MC
1297002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1298209238afSKris Buschelman 
1299209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1300209238afSKris Buschelman    and MATMPIDENSE otherwise.
1301209238afSKris Buschelman 
1302209238afSKris Buschelman    Options Database Keys:
1303209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1304209238afSKris Buschelman 
1305209238afSKris Buschelman   Level: beginner
1306209238afSKris Buschelman 
130701b82886SBarry Smith 
1308209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1309209238afSKris Buschelman M*/
1310209238afSKris Buschelman 
13114a2ae208SSatish Balay #undef __FUNCT__
13124a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1313273d9f13SBarry Smith /*@C
1314273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1315273d9f13SBarry Smith 
1316273d9f13SBarry Smith    Not collective
1317273d9f13SBarry Smith 
1318273d9f13SBarry Smith    Input Parameters:
1319273d9f13SBarry Smith .  A - the matrix
13200298fd71SBarry Smith -  data - optional location of matrix data.  Set data=NULL for PETSc
1321273d9f13SBarry Smith    to control all matrix memory allocation.
1322273d9f13SBarry Smith 
1323273d9f13SBarry Smith    Notes:
1324273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1325273d9f13SBarry Smith    storage by columns.
1326273d9f13SBarry Smith 
1327273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1328273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
13290298fd71SBarry Smith    set data=NULL.
1330273d9f13SBarry Smith 
1331273d9f13SBarry Smith    Level: intermediate
1332273d9f13SBarry Smith 
1333273d9f13SBarry Smith .keywords: matrix,dense, parallel
1334273d9f13SBarry Smith 
1335273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1336273d9f13SBarry Smith @*/
13377087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1338273d9f13SBarry Smith {
13394ac538c5SBarry Smith   PetscErrorCode ierr;
1340273d9f13SBarry Smith 
1341273d9f13SBarry Smith   PetscFunctionBegin;
13424ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(mat,data));CHKERRQ(ierr);
1343273d9f13SBarry Smith   PetscFunctionReturn(0);
1344273d9f13SBarry Smith }
1345273d9f13SBarry Smith 
13464a2ae208SSatish Balay #undef __FUNCT__
134769b1f4b7SBarry Smith #define __FUNCT__ "MatCreateDense"
13488965ea79SLois Curfman McInnes /*@C
134969b1f4b7SBarry Smith    MatCreateDense - Creates a parallel matrix in dense format.
13508965ea79SLois Curfman McInnes 
1351db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1352db81eaa0SLois Curfman McInnes 
13538965ea79SLois Curfman McInnes    Input Parameters:
1354db81eaa0SLois Curfman McInnes +  comm - MPI communicator
13558965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1356db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
13578965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1358db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
13590298fd71SBarry Smith -  data - optional location of matrix data.  Set data=NULL (NULL_SCALAR for Fortran users) for PETSc
1360dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
13618965ea79SLois Curfman McInnes 
13628965ea79SLois Curfman McInnes    Output Parameter:
1363477f1c0bSLois Curfman McInnes .  A - the matrix
13648965ea79SLois Curfman McInnes 
1365b259b22eSLois Curfman McInnes    Notes:
136639ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
136739ddd567SLois Curfman McInnes    storage by columns.
13688965ea79SLois Curfman McInnes 
136918f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
137018f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
13710298fd71SBarry Smith    set data=NULL (NULL_SCALAR for Fortran users).
137218f449edSLois Curfman McInnes 
13738965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
13748965ea79SLois Curfman McInnes    (possibly both).
13758965ea79SLois Curfman McInnes 
1376027ccd11SLois Curfman McInnes    Level: intermediate
1377027ccd11SLois Curfman McInnes 
137839ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
13798965ea79SLois Curfman McInnes 
138039ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
13818965ea79SLois Curfman McInnes @*/
138269b1f4b7SBarry Smith PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
13838965ea79SLois Curfman McInnes {
13846849ba73SBarry Smith   PetscErrorCode ierr;
138513f74950SBarry Smith   PetscMPIInt    size;
13868965ea79SLois Curfman McInnes 
13873a40ed3dSBarry Smith   PetscFunctionBegin;
1388f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1389f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1390273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1391273d9f13SBarry Smith   if (size > 1) {
1392273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1393273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1394273d9f13SBarry Smith   } else {
1395273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1396273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
13978c469469SLois Curfman McInnes   }
13983a40ed3dSBarry Smith   PetscFunctionReturn(0);
13998965ea79SLois Curfman McInnes }
14008965ea79SLois Curfman McInnes 
14014a2ae208SSatish Balay #undef __FUNCT__
14024a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
14036849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
14048965ea79SLois Curfman McInnes {
14058965ea79SLois Curfman McInnes   Mat            mat;
14063501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1407dfbe8321SBarry Smith   PetscErrorCode ierr;
14088965ea79SLois Curfman McInnes 
14093a40ed3dSBarry Smith   PetscFunctionBegin;
14108965ea79SLois Curfman McInnes   *newmat = 0;
1411*ce94432eSBarry Smith   ierr    = MatCreate(PetscObjectComm((PetscObject)A),&mat);CHKERRQ(ierr);
1412d0f46423SBarry Smith   ierr    = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
14137adad957SLisandro Dalcin   ierr    = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1414834f8fabSBarry Smith   a       = (Mat_MPIDense*)mat->data;
1415e04c1aa4SHong Zhang   ierr    = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
14165aa7edbeSHong Zhang 
1417d5f3da31SBarry Smith   mat->factortype   = A->factortype;
1418c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1419273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
14208965ea79SLois Curfman McInnes 
14218965ea79SLois Curfman McInnes   a->size         = oldmat->size;
14228965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1423e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1424b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
14253782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1426e04c1aa4SHong Zhang 
14271e1e43feSBarry Smith   ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr);
14281e1e43feSBarry Smith   ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr);
14298965ea79SLois Curfman McInnes 
1430329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
14315609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
143252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
143301b82886SBarry Smith 
14348965ea79SLois Curfman McInnes   *newmat = mat;
14353a40ed3dSBarry Smith   PetscFunctionReturn(0);
14368965ea79SLois Curfman McInnes }
14378965ea79SLois Curfman McInnes 
14384a2ae208SSatish Balay #undef __FUNCT__
14395bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
14405bba2384SShri Abhyankar PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat,PetscInt sizesset)
144186d161a7SShri Abhyankar {
144286d161a7SShri Abhyankar   PetscErrorCode ierr;
144386d161a7SShri Abhyankar   PetscMPIInt    rank,size;
144486d161a7SShri Abhyankar   PetscInt       *rowners,i,m,nz,j;
144586d161a7SShri Abhyankar   PetscScalar    *array,*vals,*vals_ptr;
144686d161a7SShri Abhyankar 
144786d161a7SShri Abhyankar   PetscFunctionBegin;
144886d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
144986d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
145086d161a7SShri Abhyankar 
145186d161a7SShri Abhyankar   /* determine ownership of all rows */
145286d161a7SShri Abhyankar   if (newmat->rmap->n < 0) m = M/size + ((M % size) > rank);
145386d161a7SShri Abhyankar   else m = newmat->rmap->n;
145486d161a7SShri Abhyankar   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
145586d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
145686d161a7SShri Abhyankar   rowners[0] = 0;
145786d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
145886d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
145986d161a7SShri Abhyankar   }
146086d161a7SShri Abhyankar 
146186d161a7SShri Abhyankar   if (!sizesset) {
146286d161a7SShri Abhyankar     ierr = MatSetSizes(newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
146386d161a7SShri Abhyankar   }
14640298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
14658c778c55SBarry Smith   ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr);
146686d161a7SShri Abhyankar 
146786d161a7SShri Abhyankar   if (!rank) {
146886d161a7SShri Abhyankar     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
146986d161a7SShri Abhyankar 
147086d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
147186d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
147286d161a7SShri Abhyankar 
147386d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
147486d161a7SShri Abhyankar     vals_ptr = vals;
147586d161a7SShri Abhyankar     for (i=0; i<m; i++) {
147686d161a7SShri Abhyankar       for (j=0; j<N; j++) {
147786d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
147886d161a7SShri Abhyankar       }
147986d161a7SShri Abhyankar     }
148086d161a7SShri Abhyankar 
148186d161a7SShri Abhyankar     /* read in other processors and ship out */
148286d161a7SShri Abhyankar     for (i=1; i<size; i++) {
148386d161a7SShri Abhyankar       nz   = (rowners[i+1] - rowners[i])*N;
148486d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1485a25532f0SBarry Smith       ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
148686d161a7SShri Abhyankar     }
148786d161a7SShri Abhyankar   } else {
148886d161a7SShri Abhyankar     /* receive numeric values */
148986d161a7SShri Abhyankar     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
149086d161a7SShri Abhyankar 
149186d161a7SShri Abhyankar     /* receive message of values*/
1492a25532f0SBarry Smith     ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
149386d161a7SShri Abhyankar 
149486d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
149586d161a7SShri Abhyankar     vals_ptr = vals;
149686d161a7SShri Abhyankar     for (i=0; i<m; i++) {
149786d161a7SShri Abhyankar       for (j=0; j<N; j++) {
149886d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
149986d161a7SShri Abhyankar       }
150086d161a7SShri Abhyankar     }
150186d161a7SShri Abhyankar   }
15028c778c55SBarry Smith   ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr);
150386d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
150486d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
150586d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150686d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150786d161a7SShri Abhyankar   PetscFunctionReturn(0);
150886d161a7SShri Abhyankar }
150986d161a7SShri Abhyankar 
151086d161a7SShri Abhyankar #undef __FUNCT__
15115bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense"
1512112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer)
151386d161a7SShri Abhyankar {
151486d161a7SShri Abhyankar   PetscScalar    *vals,*svals;
1515*ce94432eSBarry Smith   MPI_Comm       comm;
151686d161a7SShri Abhyankar   MPI_Status     status;
151786d161a7SShri Abhyankar   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
151886d161a7SShri Abhyankar   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
151986d161a7SShri Abhyankar   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
152086d161a7SShri Abhyankar   PetscInt       i,nz,j,rstart,rend,sizesset=1,grows,gcols;
152186d161a7SShri Abhyankar   int            fd;
152286d161a7SShri Abhyankar   PetscErrorCode ierr;
152386d161a7SShri Abhyankar 
152486d161a7SShri Abhyankar   PetscFunctionBegin;
1525*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
152686d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
152786d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
152886d161a7SShri Abhyankar   if (!rank) {
152986d161a7SShri Abhyankar     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
153086d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,(char*)header,4,PETSC_INT);CHKERRQ(ierr);
153186d161a7SShri Abhyankar     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
153286d161a7SShri Abhyankar   }
153386d161a7SShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) sizesset = 0;
153486d161a7SShri Abhyankar 
153586d161a7SShri Abhyankar   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
153686d161a7SShri Abhyankar   M    = header[1]; N = header[2]; nz = header[3];
153786d161a7SShri Abhyankar 
153886d161a7SShri Abhyankar   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
153986d161a7SShri Abhyankar   if (sizesset && newmat->rmap->N < 0) newmat->rmap->N = M;
154086d161a7SShri Abhyankar   if (sizesset && newmat->cmap->N < 0) newmat->cmap->N = N;
154186d161a7SShri Abhyankar 
154286d161a7SShri Abhyankar   /* If global sizes are set, check if they are consistent with that given in the file */
154386d161a7SShri Abhyankar   if (sizesset) {
154486d161a7SShri Abhyankar     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
154586d161a7SShri Abhyankar   }
1546abd38a8fSBarry Smith   if (sizesset && newmat->rmap->N != grows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%d) and input matrix has (%d)",M,grows);
1547abd38a8fSBarry Smith   if (sizesset && newmat->cmap->N != gcols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%d) and input matrix has (%d)",N,gcols);
154886d161a7SShri Abhyankar 
154986d161a7SShri Abhyankar   /*
155086d161a7SShri Abhyankar        Handle case where matrix is stored on disk as a dense matrix
155186d161a7SShri Abhyankar   */
155286d161a7SShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
15535bba2384SShri Abhyankar     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat,sizesset);CHKERRQ(ierr);
155486d161a7SShri Abhyankar     PetscFunctionReturn(0);
155586d161a7SShri Abhyankar   }
155686d161a7SShri Abhyankar 
155786d161a7SShri Abhyankar   /* determine ownership of all rows */
15582205254eSKarl Rupp   if (newmat->rmap->n < 0) {
15592205254eSKarl Rupp     ierr = PetscMPIIntCast(M/size + ((M % size) > rank),&m);CHKERRQ(ierr);
15602205254eSKarl Rupp   } else {
15612205254eSKarl Rupp     ierr = PetscMPIIntCast(newmat->rmap->n,&m);CHKERRQ(ierr);
15622205254eSKarl Rupp   }
156386d161a7SShri Abhyankar   ierr       = PetscMalloc((size+2)*sizeof(PetscMPIInt),&rowners);CHKERRQ(ierr);
156486d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
156586d161a7SShri Abhyankar   rowners[0] = 0;
156686d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
156786d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
156886d161a7SShri Abhyankar   }
156986d161a7SShri Abhyankar   rstart = rowners[rank];
157086d161a7SShri Abhyankar   rend   = rowners[rank+1];
157186d161a7SShri Abhyankar 
157286d161a7SShri Abhyankar   /* distribute row lengths to all processors */
157386d161a7SShri Abhyankar   ierr = PetscMalloc2(rend-rstart,PetscInt,&ourlens,rend-rstart,PetscInt,&offlens);CHKERRQ(ierr);
157486d161a7SShri Abhyankar   if (!rank) {
157586d161a7SShri Abhyankar     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
157686d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
157786d161a7SShri Abhyankar     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
157886d161a7SShri Abhyankar     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
157986d161a7SShri Abhyankar     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
158086d161a7SShri Abhyankar     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
158186d161a7SShri Abhyankar   } else {
158286d161a7SShri Abhyankar     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
158386d161a7SShri Abhyankar   }
158486d161a7SShri Abhyankar 
158586d161a7SShri Abhyankar   if (!rank) {
158686d161a7SShri Abhyankar     /* calculate the number of nonzeros on each processor */
158786d161a7SShri Abhyankar     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
158886d161a7SShri Abhyankar     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
158986d161a7SShri Abhyankar     for (i=0; i<size; i++) {
159086d161a7SShri Abhyankar       for (j=rowners[i]; j< rowners[i+1]; j++) {
159186d161a7SShri Abhyankar         procsnz[i] += rowlengths[j];
159286d161a7SShri Abhyankar       }
159386d161a7SShri Abhyankar     }
159486d161a7SShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
159586d161a7SShri Abhyankar 
159686d161a7SShri Abhyankar     /* determine max buffer needed and allocate it */
159786d161a7SShri Abhyankar     maxnz = 0;
159886d161a7SShri Abhyankar     for (i=0; i<size; i++) {
159986d161a7SShri Abhyankar       maxnz = PetscMax(maxnz,procsnz[i]);
160086d161a7SShri Abhyankar     }
160186d161a7SShri Abhyankar     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
160286d161a7SShri Abhyankar 
160386d161a7SShri Abhyankar     /* read in my part of the matrix column indices  */
160486d161a7SShri Abhyankar     nz   = procsnz[0];
160586d161a7SShri Abhyankar     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
160686d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
160786d161a7SShri Abhyankar 
160886d161a7SShri Abhyankar     /* read in every one elses and ship off */
160986d161a7SShri Abhyankar     for (i=1; i<size; i++) {
161086d161a7SShri Abhyankar       nz   = procsnz[i];
161186d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
161286d161a7SShri Abhyankar       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
161386d161a7SShri Abhyankar     }
161486d161a7SShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
161586d161a7SShri Abhyankar   } else {
161686d161a7SShri Abhyankar     /* determine buffer space needed for message */
161786d161a7SShri Abhyankar     nz = 0;
161886d161a7SShri Abhyankar     for (i=0; i<m; i++) {
161986d161a7SShri Abhyankar       nz += ourlens[i];
162086d161a7SShri Abhyankar     }
162186d161a7SShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
162286d161a7SShri Abhyankar 
162386d161a7SShri Abhyankar     /* receive message of column indices*/
162486d161a7SShri Abhyankar     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
162586d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
162686d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
162786d161a7SShri Abhyankar   }
162886d161a7SShri Abhyankar 
162986d161a7SShri Abhyankar   /* loop over local rows, determining number of off diagonal entries */
163086d161a7SShri Abhyankar   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
163186d161a7SShri Abhyankar   jj   = 0;
163286d161a7SShri Abhyankar   for (i=0; i<m; i++) {
163386d161a7SShri Abhyankar     for (j=0; j<ourlens[i]; j++) {
163486d161a7SShri Abhyankar       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
163586d161a7SShri Abhyankar       jj++;
163686d161a7SShri Abhyankar     }
163786d161a7SShri Abhyankar   }
163886d161a7SShri Abhyankar 
163986d161a7SShri Abhyankar   /* create our matrix */
16402205254eSKarl Rupp   for (i=0; i<m; i++) ourlens[i] -= offlens[i];
164186d161a7SShri Abhyankar 
164286d161a7SShri Abhyankar   if (!sizesset) {
164386d161a7SShri Abhyankar     ierr = MatSetSizes(newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
164486d161a7SShri Abhyankar   }
16450298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
16462205254eSKarl Rupp   for (i=0; i<m; i++) ourlens[i] += offlens[i];
164786d161a7SShri Abhyankar 
164886d161a7SShri Abhyankar   if (!rank) {
164986d161a7SShri Abhyankar     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
165086d161a7SShri Abhyankar 
165186d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
165286d161a7SShri Abhyankar     nz   = procsnz[0];
165386d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
165486d161a7SShri Abhyankar 
165586d161a7SShri Abhyankar     /* insert into matrix */
165686d161a7SShri Abhyankar     jj      = rstart;
165786d161a7SShri Abhyankar     smycols = mycols;
165886d161a7SShri Abhyankar     svals   = vals;
165986d161a7SShri Abhyankar     for (i=0; i<m; i++) {
166086d161a7SShri Abhyankar       ierr     = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
166186d161a7SShri Abhyankar       smycols += ourlens[i];
166286d161a7SShri Abhyankar       svals   += ourlens[i];
166386d161a7SShri Abhyankar       jj++;
166486d161a7SShri Abhyankar     }
166586d161a7SShri Abhyankar 
166686d161a7SShri Abhyankar     /* read in other processors and ship out */
166786d161a7SShri Abhyankar     for (i=1; i<size; i++) {
166886d161a7SShri Abhyankar       nz   = procsnz[i];
166986d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
167086d161a7SShri Abhyankar       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr);
167186d161a7SShri Abhyankar     }
167286d161a7SShri Abhyankar     ierr = PetscFree(procsnz);CHKERRQ(ierr);
167386d161a7SShri Abhyankar   } else {
167486d161a7SShri Abhyankar     /* receive numeric values */
167586d161a7SShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
167686d161a7SShri Abhyankar 
167786d161a7SShri Abhyankar     /* receive message of values*/
167886d161a7SShri Abhyankar     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr);
167986d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
168086d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
168186d161a7SShri Abhyankar 
168286d161a7SShri Abhyankar     /* insert into matrix */
168386d161a7SShri Abhyankar     jj      = rstart;
168486d161a7SShri Abhyankar     smycols = mycols;
168586d161a7SShri Abhyankar     svals   = vals;
168686d161a7SShri Abhyankar     for (i=0; i<m; i++) {
168786d161a7SShri Abhyankar       ierr     = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
168886d161a7SShri Abhyankar       smycols += ourlens[i];
168986d161a7SShri Abhyankar       svals   += ourlens[i];
169086d161a7SShri Abhyankar       jj++;
169186d161a7SShri Abhyankar     }
169286d161a7SShri Abhyankar   }
169386d161a7SShri Abhyankar   ierr = PetscFree2(ourlens,offlens);CHKERRQ(ierr);
169486d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
169586d161a7SShri Abhyankar   ierr = PetscFree(mycols);CHKERRQ(ierr);
169686d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
169786d161a7SShri Abhyankar 
169886d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169986d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
170086d161a7SShri Abhyankar   PetscFunctionReturn(0);
170186d161a7SShri Abhyankar }
170286d161a7SShri Abhyankar 
170386d161a7SShri Abhyankar #undef __FUNCT__
17046e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense"
1705ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
17066e4ee0c6SHong Zhang {
17076e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
17086e4ee0c6SHong Zhang   Mat            a,b;
1709ace3abfcSBarry Smith   PetscBool      flg;
17106e4ee0c6SHong Zhang   PetscErrorCode ierr;
171190ace30eSBarry Smith 
17126e4ee0c6SHong Zhang   PetscFunctionBegin;
17136e4ee0c6SHong Zhang   a    = matA->A;
17146e4ee0c6SHong Zhang   b    = matB->A;
17156e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
1716*ce94432eSBarry Smith   ierr = MPI_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
17176e4ee0c6SHong Zhang   PetscFunctionReturn(0);
17186e4ee0c6SHong Zhang }
171990ace30eSBarry Smith 
1720