xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 8572280aa1ec174234a4e348c599e5742a8e2c3f)
1be1d678aSKris Buschelman 
2ed3cc1f0SBarry Smith /*
3ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
4ed3cc1f0SBarry Smith */
5ed3cc1f0SBarry Smith 
604fea9ffSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>    /*I   "petscmat.h"  I*/
88949adfdSHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>
9baa3c1c6SHong Zhang #include <petscblaslapack.h>
108965ea79SLois Curfman McInnes 
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 
38ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
39ba8c8a56SBarry Smith {
40ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
41ba8c8a56SBarry Smith   PetscErrorCode ierr;
42d0f46423SBarry Smith   PetscInt       lrow,rstart = A->rmap->rstart,rend = A->rmap->rend;
43ba8c8a56SBarry Smith 
44ba8c8a56SBarry Smith   PetscFunctionBegin;
45e7e72b3dSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows");
46ba8c8a56SBarry Smith   lrow = row - rstart;
47ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt**)idx,(const PetscScalar**)v);CHKERRQ(ierr);
48ba8c8a56SBarry Smith   PetscFunctionReturn(0);
49ba8c8a56SBarry Smith }
50ba8c8a56SBarry Smith 
51ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
52ba8c8a56SBarry Smith {
53ba8c8a56SBarry Smith   PetscErrorCode ierr;
54ba8c8a56SBarry Smith 
55ba8c8a56SBarry Smith   PetscFunctionBegin;
56ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
57ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
58ba8c8a56SBarry Smith   PetscFunctionReturn(0);
59ba8c8a56SBarry Smith }
60ba8c8a56SBarry Smith 
6111bd1e4dSLisandro Dalcin PetscErrorCode  MatGetDiagonalBlock_MPIDense(Mat A,Mat *a)
620de54da6SSatish Balay {
630de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
646849ba73SBarry Smith   PetscErrorCode ierr;
65d0f46423SBarry Smith   PetscInt       m = A->rmap->n,rstart = A->rmap->rstart;
6687828ca2SBarry Smith   PetscScalar    *array;
670de54da6SSatish Balay   MPI_Comm       comm;
6811bd1e4dSLisandro Dalcin   Mat            B;
690de54da6SSatish Balay 
700de54da6SSatish Balay   PetscFunctionBegin;
71e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only square matrices supported.");
720de54da6SSatish Balay 
7311bd1e4dSLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);CHKERRQ(ierr);
7411bd1e4dSLisandro Dalcin   if (!B) {
750de54da6SSatish Balay     ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
7611bd1e4dSLisandro Dalcin     ierr = MatCreate(comm,&B);CHKERRQ(ierr);
7711bd1e4dSLisandro Dalcin     ierr = MatSetSizes(B,m,m,m,m);CHKERRQ(ierr);
7811bd1e4dSLisandro Dalcin     ierr = MatSetType(B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr);
798c778c55SBarry Smith     ierr = MatDenseGetArray(mdn->A,&array);CHKERRQ(ierr);
8011bd1e4dSLisandro Dalcin     ierr = MatSeqDenseSetPreallocation(B,array+m*rstart);CHKERRQ(ierr);
818c778c55SBarry Smith     ierr = MatDenseRestoreArray(mdn->A,&array);CHKERRQ(ierr);
8211bd1e4dSLisandro Dalcin     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8311bd1e4dSLisandro Dalcin     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8411bd1e4dSLisandro Dalcin     ierr = PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);CHKERRQ(ierr);
8511bd1e4dSLisandro Dalcin     *a   = B;
8611bd1e4dSLisandro Dalcin     ierr = MatDestroy(&B);CHKERRQ(ierr);
872205254eSKarl Rupp   } else *a = B;
880de54da6SSatish Balay   PetscFunctionReturn(0);
890de54da6SSatish Balay }
900de54da6SSatish Balay 
9113f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
928965ea79SLois Curfman McInnes {
9339b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
94dfbe8321SBarry Smith   PetscErrorCode ierr;
95d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
96ace3abfcSBarry Smith   PetscBool      roworiented = A->roworiented;
978965ea79SLois Curfman McInnes 
983a40ed3dSBarry Smith   PetscFunctionBegin;
998965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1005ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
101e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1028965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1038965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
10439b7565bSBarry Smith       if (roworiented) {
10539b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1063a40ed3dSBarry Smith       } else {
1078965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1085ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
109e32f2f54SBarry Smith           if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
11039b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
11139b7565bSBarry Smith         }
1128965ea79SLois Curfman McInnes       }
1132205254eSKarl Rupp     } else if (!A->donotstash) {
1145080c13bSMatthew G Knepley       mat->assembled = PETSC_FALSE;
11539b7565bSBarry Smith       if (roworiented) {
116b400d20cSBarry Smith         ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n,PETSC_FALSE);CHKERRQ(ierr);
117d36fbae8SSatish Balay       } else {
118b400d20cSBarry Smith         ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m,PETSC_FALSE);CHKERRQ(ierr);
11939b7565bSBarry Smith       }
120b49de8d1SLois Curfman McInnes     }
121b49de8d1SLois Curfman McInnes   }
1223a40ed3dSBarry Smith   PetscFunctionReturn(0);
123b49de8d1SLois Curfman McInnes }
124b49de8d1SLois Curfman McInnes 
12513f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
126b49de8d1SLois Curfman McInnes {
127b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
128dfbe8321SBarry Smith   PetscErrorCode ierr;
129d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
130b49de8d1SLois Curfman McInnes 
1313a40ed3dSBarry Smith   PetscFunctionBegin;
132b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
133e32f2f54SBarry Smith     if (idxm[i] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
134e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
135b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
136b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
137b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
138e32f2f54SBarry Smith         if (idxn[j] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
139e7e72b3dSBarry Smith         if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
140b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
141b49de8d1SLois Curfman McInnes       }
142e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
1438965ea79SLois Curfman McInnes   }
1443a40ed3dSBarry Smith   PetscFunctionReturn(0);
1458965ea79SLois Curfman McInnes }
1468965ea79SLois Curfman McInnes 
147d3042a70SBarry Smith static PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[])
148ff14e315SSatish Balay {
149ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
150dfbe8321SBarry Smith   PetscErrorCode ierr;
151ff14e315SSatish Balay 
1523a40ed3dSBarry Smith   PetscFunctionBegin;
1538c778c55SBarry Smith   ierr = MatDenseGetArray(a->A,array);CHKERRQ(ierr);
1543a40ed3dSBarry Smith   PetscFunctionReturn(0);
155ff14e315SSatish Balay }
156ff14e315SSatish Balay 
157*8572280aSBarry Smith static PetscErrorCode MatDenseGetArrayRead_MPIDense(Mat A,const PetscScalar *array[])
158*8572280aSBarry Smith {
159*8572280aSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
160*8572280aSBarry Smith   PetscErrorCode ierr;
161*8572280aSBarry Smith 
162*8572280aSBarry Smith   PetscFunctionBegin;
163*8572280aSBarry Smith   ierr = MatDenseGetArrayRead(a->A,array);CHKERRQ(ierr);
164*8572280aSBarry Smith   PetscFunctionReturn(0);
165*8572280aSBarry Smith }
166*8572280aSBarry Smith 
167d3042a70SBarry Smith static PetscErrorCode MatDensePlaceArray_MPIDense(Mat A,const PetscScalar array[])
168d3042a70SBarry Smith {
169d3042a70SBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
170d3042a70SBarry Smith   PetscErrorCode ierr;
171d3042a70SBarry Smith 
172d3042a70SBarry Smith   PetscFunctionBegin;
173d3042a70SBarry Smith   ierr = MatDensePlaceArray(a->A,array);CHKERRQ(ierr);
174d3042a70SBarry Smith   PetscFunctionReturn(0);
175d3042a70SBarry Smith }
176d3042a70SBarry Smith 
177d3042a70SBarry Smith static PetscErrorCode MatDenseResetArray_MPIDense(Mat A)
178d3042a70SBarry Smith {
179d3042a70SBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
180d3042a70SBarry Smith   PetscErrorCode ierr;
181d3042a70SBarry Smith 
182d3042a70SBarry Smith   PetscFunctionBegin;
183d3042a70SBarry Smith   ierr = MatDenseResetArray(a->A);CHKERRQ(ierr);
184d3042a70SBarry Smith   PetscFunctionReturn(0);
185d3042a70SBarry Smith }
186d3042a70SBarry Smith 
1877dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
188ca3fa75bSLois Curfman McInnes {
189ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat  = (Mat_MPIDense*)A->data,*newmatd;
190ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
1916849ba73SBarry Smith   PetscErrorCode ierr;
1924aa3045dSJed Brown   PetscInt       i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols;
1935d0c19d7SBarry Smith   const PetscInt *irow,*icol;
19487828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
195ca3fa75bSLois Curfman McInnes   Mat            newmat;
1964aa3045dSJed Brown   IS             iscol_local;
197ca3fa75bSLois Curfman McInnes 
198ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
1994aa3045dSJed Brown   ierr = ISAllGather(iscol,&iscol_local);CHKERRQ(ierr);
200ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
2014aa3045dSJed Brown   ierr = ISGetIndices(iscol_local,&icol);CHKERRQ(ierr);
202b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
203b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
2044aa3045dSJed Brown   ierr = ISGetSize(iscol,&Ncols);CHKERRQ(ierr); /* global number of columns, size of iscol_local */
205ca3fa75bSLois Curfman McInnes 
206ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
2077eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
2087eba5e9cSLois Curfman McInnes      original matrix! */
209ca3fa75bSLois Curfman McInnes 
210ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
2117eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
212ca3fa75bSLois Curfman McInnes 
213ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
214ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
215e32f2f54SBarry Smith     /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2167eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
217ca3fa75bSLois Curfman McInnes     newmat = *B;
218ca3fa75bSLois Curfman McInnes   } else {
219ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
220ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr);
2214aa3045dSJed Brown     ierr = MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);CHKERRQ(ierr);
2227adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
2230298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
224ca3fa75bSLois Curfman McInnes   }
225ca3fa75bSLois Curfman McInnes 
226ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
227ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
228ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense*)newmatd->A->data)->v;
229ca3fa75bSLois Curfman McInnes 
2304aa3045dSJed Brown   for (i=0; i<Ncols; i++) {
23125a33276SHong Zhang     av = v + ((Mat_SeqDense*)mat->A->data)->lda*icol[i];
232ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2337eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
234ca3fa75bSLois Curfman McInnes     }
235ca3fa75bSLois Curfman McInnes   }
236ca3fa75bSLois Curfman McInnes 
237ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
238ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
239ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
240ca3fa75bSLois Curfman McInnes 
241ca3fa75bSLois Curfman McInnes   /* Free work space */
242ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
2435bdf786aSShri Abhyankar   ierr = ISRestoreIndices(iscol_local,&icol);CHKERRQ(ierr);
24432bb1f2dSLisandro Dalcin   ierr = ISDestroy(&iscol_local);CHKERRQ(ierr);
245ca3fa75bSLois Curfman McInnes   *B   = newmat;
246ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
247ca3fa75bSLois Curfman McInnes }
248ca3fa75bSLois Curfman McInnes 
2498c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[])
250ff14e315SSatish Balay {
25173a71a0fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
25273a71a0fSBarry Smith   PetscErrorCode ierr;
25373a71a0fSBarry Smith 
2543a40ed3dSBarry Smith   PetscFunctionBegin;
2558c778c55SBarry Smith   ierr = MatDenseRestoreArray(a->A,array);CHKERRQ(ierr);
2563a40ed3dSBarry Smith   PetscFunctionReturn(0);
257ff14e315SSatish Balay }
258ff14e315SSatish Balay 
259*8572280aSBarry Smith PetscErrorCode MatDenseRestoreArrayRead_MPIDense(Mat A,const PetscScalar *array[])
260*8572280aSBarry Smith {
261*8572280aSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
262*8572280aSBarry Smith   PetscErrorCode ierr;
263*8572280aSBarry Smith 
264*8572280aSBarry Smith   PetscFunctionBegin;
265*8572280aSBarry Smith   ierr = MatDenseRestoreArrayRead(a->A,array);CHKERRQ(ierr);
266*8572280aSBarry Smith   PetscFunctionReturn(0);
267*8572280aSBarry Smith }
268*8572280aSBarry Smith 
269dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2708965ea79SLois Curfman McInnes {
27139ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
272ce94432eSBarry Smith   MPI_Comm       comm;
273dfbe8321SBarry Smith   PetscErrorCode ierr;
27413f74950SBarry Smith   PetscInt       nstash,reallocs;
2758965ea79SLois Curfman McInnes   InsertMode     addv;
2768965ea79SLois Curfman McInnes 
2773a40ed3dSBarry Smith   PetscFunctionBegin;
278ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
2798965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
280b2566f29SBarry Smith   ierr = MPIU_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);CHKERRQ(ierr);
281ce94432eSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
282e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2838965ea79SLois Curfman McInnes 
284d0f46423SBarry Smith   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr);
2858798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
286ae15b995SBarry Smith   ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
2873a40ed3dSBarry Smith   PetscFunctionReturn(0);
2888965ea79SLois Curfman McInnes }
2898965ea79SLois Curfman McInnes 
290dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2918965ea79SLois Curfman McInnes {
29239ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn=(Mat_MPIDense*)mat->data;
2936849ba73SBarry Smith   PetscErrorCode ierr;
29413f74950SBarry Smith   PetscInt       i,*row,*col,flg,j,rstart,ncols;
29513f74950SBarry Smith   PetscMPIInt    n;
29687828ca2SBarry Smith   PetscScalar    *val;
2978965ea79SLois Curfman McInnes 
2983a40ed3dSBarry Smith   PetscFunctionBegin;
2998965ea79SLois Curfman McInnes   /*  wait on receives */
3007ef1d9bdSSatish Balay   while (1) {
3018798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
3027ef1d9bdSSatish Balay     if (!flg) break;
3038965ea79SLois Curfman McInnes 
3047ef1d9bdSSatish Balay     for (i=0; i<n;) {
3057ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
3062205254eSKarl Rupp       for (j=i,rstart=row[j]; j<n; j++) {
3072205254eSKarl Rupp         if (row[j] != rstart) break;
3082205254eSKarl Rupp       }
3097ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
3107ef1d9bdSSatish Balay       else       ncols = n-i;
3117ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
3124b4eb8d3SJed Brown       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);CHKERRQ(ierr);
3137ef1d9bdSSatish Balay       i    = j;
3148965ea79SLois Curfman McInnes     }
3157ef1d9bdSSatish Balay   }
3168798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
3178965ea79SLois Curfman McInnes 
31839ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
31939ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3208965ea79SLois Curfman McInnes 
3216d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
32239ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3238965ea79SLois Curfman McInnes   }
3243a40ed3dSBarry Smith   PetscFunctionReturn(0);
3258965ea79SLois Curfman McInnes }
3268965ea79SLois Curfman McInnes 
327dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3288965ea79SLois Curfman McInnes {
329dfbe8321SBarry Smith   PetscErrorCode ierr;
33039ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3313a40ed3dSBarry Smith 
3323a40ed3dSBarry Smith   PetscFunctionBegin;
3333a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3343a40ed3dSBarry Smith   PetscFunctionReturn(0);
3358965ea79SLois Curfman McInnes }
3368965ea79SLois Curfman McInnes 
3378965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3388965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3398965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3403501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3418965ea79SLois Curfman McInnes    routine.
3428965ea79SLois Curfman McInnes */
3432b40b63fSBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
3448965ea79SLois Curfman McInnes {
34539ddd567SLois Curfman McInnes   Mat_MPIDense      *l = (Mat_MPIDense*)A->data;
3466849ba73SBarry Smith   PetscErrorCode    ierr;
347d0f46423SBarry Smith   PetscInt          i,*owners = A->rmap->range;
34876ec1555SBarry Smith   PetscInt          *sizes,j,idx,nsends;
34913f74950SBarry Smith   PetscInt          nmax,*svalues,*starts,*owner,nrecvs;
3507adad957SLisandro Dalcin   PetscInt          *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source;
35113f74950SBarry Smith   PetscInt          *lens,*lrows,*values;
35213f74950SBarry Smith   PetscMPIInt       n,imdex,rank = l->rank,size = l->size;
353ce94432eSBarry Smith   MPI_Comm          comm;
3548965ea79SLois Curfman McInnes   MPI_Request       *send_waits,*recv_waits;
3558965ea79SLois Curfman McInnes   MPI_Status        recv_status,*send_status;
356ace3abfcSBarry Smith   PetscBool         found;
35797b48c8fSBarry Smith   const PetscScalar *xx;
35897b48c8fSBarry Smith   PetscScalar       *bb;
3598965ea79SLois Curfman McInnes 
3603a40ed3dSBarry Smith   PetscFunctionBegin;
361ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
362ce94432eSBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(comm,PETSC_ERR_SUP,"Only handles square matrices");
363b9679d65SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership");
3648965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
365f628708eSJed Brown   ierr = PetscCalloc1(2*size,&sizes);CHKERRQ(ierr);
366037dbc42SBarry Smith   ierr = PetscMalloc1(N+1,&owner);CHKERRQ(ierr);  /* see note*/
3678965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3688965ea79SLois Curfman McInnes     idx   = rows[i];
36935d8aa7fSBarry Smith     found = PETSC_FALSE;
3708965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3718965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
37276ec1555SBarry Smith         sizes[2*j]++; sizes[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3738965ea79SLois Curfman McInnes       }
3748965ea79SLois Curfman McInnes     }
375e32f2f54SBarry Smith     if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3768965ea79SLois Curfman McInnes   }
3772205254eSKarl Rupp   nsends = 0;
37876ec1555SBarry Smith   for (i=0; i<size; i++) nsends += sizes[2*i+1];
3798965ea79SLois Curfman McInnes 
3808965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
38176ec1555SBarry Smith   ierr = PetscMaxSum(comm,sizes,&nmax,&nrecvs);CHKERRQ(ierr);
3828965ea79SLois Curfman McInnes 
3838965ea79SLois Curfman McInnes   /* post receives:   */
384785e854fSJed Brown   ierr = PetscMalloc1((nrecvs+1)*(nmax+1),&rvalues);CHKERRQ(ierr);
385854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
3868965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
38713f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3888965ea79SLois Curfman McInnes   }
3898965ea79SLois Curfman McInnes 
3908965ea79SLois Curfman McInnes   /* do sends:
3918965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3928965ea79SLois Curfman McInnes          the ith processor
3938965ea79SLois Curfman McInnes   */
394854ce69bSBarry Smith   ierr = PetscMalloc1(N+1,&svalues);CHKERRQ(ierr);
395854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
396854ce69bSBarry Smith   ierr = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
3978965ea79SLois Curfman McInnes 
3988965ea79SLois Curfman McInnes   starts[0] = 0;
39976ec1555SBarry Smith   for (i=1; i<size; i++) starts[i] = starts[i-1] + sizes[2*i-2];
4002205254eSKarl Rupp   for (i=0; i<N; i++) svalues[starts[owner[i]]++] = rows[i];
4012205254eSKarl Rupp 
4022205254eSKarl Rupp   starts[0] = 0;
40376ec1555SBarry Smith   for (i=1; i<size+1; i++) starts[i] = starts[i-1] + sizes[2*i-2];
4048965ea79SLois Curfman McInnes   count = 0;
4058965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
40676ec1555SBarry Smith     if (sizes[2*i+1]) {
40776ec1555SBarry Smith       ierr = MPI_Isend(svalues+starts[i],sizes[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
4088965ea79SLois Curfman McInnes     }
4098965ea79SLois Curfman McInnes   }
410606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
4118965ea79SLois Curfman McInnes 
4128965ea79SLois Curfman McInnes   base = owners[rank];
4138965ea79SLois Curfman McInnes 
4148965ea79SLois Curfman McInnes   /*  wait on receives */
415dcca6d9dSJed Brown   ierr  = PetscMalloc2(nrecvs,&lens,nrecvs,&source);CHKERRQ(ierr);
41674ed9c26SBarry Smith   count = nrecvs;
41774ed9c26SBarry Smith   slen  = 0;
4188965ea79SLois Curfman McInnes   while (count) {
419ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4208965ea79SLois Curfman McInnes     /* unpack receives into our local space */
42113f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4222205254eSKarl Rupp 
4238965ea79SLois Curfman McInnes     source[imdex] = recv_status.MPI_SOURCE;
4248965ea79SLois Curfman McInnes     lens[imdex]   = n;
4258965ea79SLois Curfman McInnes     slen += n;
4268965ea79SLois Curfman McInnes     count--;
4278965ea79SLois Curfman McInnes   }
428606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4298965ea79SLois Curfman McInnes 
4308965ea79SLois Curfman McInnes   /* move the data into the send scatter */
431854ce69bSBarry Smith   ierr  = PetscMalloc1(slen+1,&lrows);CHKERRQ(ierr);
4328965ea79SLois Curfman McInnes   count = 0;
4338965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4348965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4358965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4368965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4378965ea79SLois Curfman McInnes     }
4388965ea79SLois Curfman McInnes   }
439606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
44074ed9c26SBarry Smith   ierr = PetscFree2(lens,source);CHKERRQ(ierr);
441606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
44276ec1555SBarry Smith   ierr = PetscFree(sizes);CHKERRQ(ierr);
4438965ea79SLois Curfman McInnes 
44497b48c8fSBarry Smith   /* fix right hand side if needed */
44597b48c8fSBarry Smith   if (x && b) {
44697b48c8fSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
44797b48c8fSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
44897b48c8fSBarry Smith     for (i=0; i<slen; i++) {
44997b48c8fSBarry Smith       bb[lrows[i]] = diag*xx[lrows[i]];
45097b48c8fSBarry Smith     }
45197b48c8fSBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
45297b48c8fSBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
45397b48c8fSBarry Smith   }
45497b48c8fSBarry Smith 
4558965ea79SLois Curfman McInnes   /* actually zap the local rows */
456b9679d65SBarry Smith   ierr = MatZeroRows(l->A,slen,lrows,0.0,0,0);CHKERRQ(ierr);
457e2eb51b1SBarry Smith   if (diag != 0.0) {
458b9679d65SBarry Smith     Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data;
459b9679d65SBarry Smith     PetscInt     m   = ll->lda, i;
460b9679d65SBarry Smith 
461b9679d65SBarry Smith     for (i=0; i<slen; i++) {
462b9679d65SBarry Smith       ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag;
463b9679d65SBarry Smith     }
464b9679d65SBarry Smith   }
465606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4668965ea79SLois Curfman McInnes 
4678965ea79SLois Curfman McInnes   /* wait on sends */
4688965ea79SLois Curfman McInnes   if (nsends) {
469785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
470ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
471606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4728965ea79SLois Curfman McInnes   }
473606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
474606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4753a40ed3dSBarry Smith   PetscFunctionReturn(0);
4768965ea79SLois Curfman McInnes }
4778965ea79SLois Curfman McInnes 
478cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
479cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
480cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
481cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);
482cc2e6a90SBarry Smith 
483dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
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);
491cc2e6a90SBarry Smith   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4923a40ed3dSBarry Smith   PetscFunctionReturn(0);
4938965ea79SLois Curfman McInnes }
4948965ea79SLois Curfman McInnes 
495dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4968965ea79SLois Curfman McInnes {
49739ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
498dfbe8321SBarry Smith   PetscErrorCode ierr;
499c456f294SBarry Smith 
5003a40ed3dSBarry Smith   PetscFunctionBegin;
501ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
502ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
503cc2e6a90SBarry Smith   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
5043a40ed3dSBarry Smith   PetscFunctionReturn(0);
5058965ea79SLois Curfman McInnes }
5068965ea79SLois Curfman McInnes 
507dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
508096963f5SLois Curfman McInnes {
509096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
510dfbe8321SBarry Smith   PetscErrorCode ierr;
51187828ca2SBarry Smith   PetscScalar    zero = 0.0;
512096963f5SLois Curfman McInnes 
5133a40ed3dSBarry Smith   PetscFunctionBegin;
5142dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
515cc2e6a90SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
516ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
517ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5183a40ed3dSBarry Smith   PetscFunctionReturn(0);
519096963f5SLois Curfman McInnes }
520096963f5SLois Curfman McInnes 
521dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
522096963f5SLois Curfman McInnes {
523096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
524dfbe8321SBarry Smith   PetscErrorCode ierr;
525096963f5SLois Curfman McInnes 
5263a40ed3dSBarry Smith   PetscFunctionBegin;
5273501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
528cc2e6a90SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
529ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
530ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5313a40ed3dSBarry Smith   PetscFunctionReturn(0);
532096963f5SLois Curfman McInnes }
533096963f5SLois Curfman McInnes 
534dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5358965ea79SLois Curfman McInnes {
53639ddd567SLois Curfman McInnes   Mat_MPIDense   *a    = (Mat_MPIDense*)A->data;
537096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
538dfbe8321SBarry Smith   PetscErrorCode ierr;
539d0f46423SBarry Smith   PetscInt       len,i,n,m = A->rmap->n,radd;
54087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
541ed3cc1f0SBarry Smith 
5423a40ed3dSBarry Smith   PetscFunctionBegin;
5432dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5441ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
545096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
546e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
547d0f46423SBarry Smith   len  = PetscMin(a->A->rmap->n,a->A->cmap->n);
548d0f46423SBarry Smith   radd = A->rmap->rstart*m;
54944cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
550096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
551096963f5SLois Curfman McInnes   }
5521ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5533a40ed3dSBarry Smith   PetscFunctionReturn(0);
5548965ea79SLois Curfman McInnes }
5558965ea79SLois Curfman McInnes 
556dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5578965ea79SLois Curfman McInnes {
5583501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
559dfbe8321SBarry Smith   PetscErrorCode ierr;
560ed3cc1f0SBarry Smith 
5613a40ed3dSBarry Smith   PetscFunctionBegin;
562aa482453SBarry Smith #if defined(PETSC_USE_LOG)
563d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N);
5648965ea79SLois Curfman McInnes #endif
5658798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
5666bf464f9SBarry Smith   ierr = MatDestroy(&mdn->A);CHKERRQ(ierr);
5676bf464f9SBarry Smith   ierr = VecDestroy(&mdn->lvec);CHKERRQ(ierr);
5686bf464f9SBarry Smith   ierr = VecScatterDestroy(&mdn->Mvctx);CHKERRQ(ierr);
56901b82886SBarry Smith 
570bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
571dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
5728baccfbdSHong Zhang 
5738baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);CHKERRQ(ierr);
574*8572280aSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);CHKERRQ(ierr);
575*8572280aSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",NULL);CHKERRQ(ierr);
576*8572280aSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",NULL);CHKERRQ(ierr);
577d3042a70SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);CHKERRQ(ierr);
578d3042a70SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);CHKERRQ(ierr);
5798baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
5808baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",NULL);CHKERRQ(ierr);
5818baccfbdSHong Zhang #endif
582bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",NULL);CHKERRQ(ierr);
583bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
584bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
585bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5868baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5878baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5888baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr);
5893a40ed3dSBarry Smith   PetscFunctionReturn(0);
5908965ea79SLois Curfman McInnes }
59139ddd567SLois Curfman McInnes 
5926849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5938965ea79SLois Curfman McInnes {
59439ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
595dfbe8321SBarry Smith   PetscErrorCode    ierr;
596aa05aa95SBarry Smith   PetscViewerFormat format;
597aa05aa95SBarry Smith   int               fd;
598d0f46423SBarry Smith   PetscInt          header[4],mmax,N = mat->cmap->N,i,j,m,k;
599aa05aa95SBarry Smith   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
600578230a0SSatish Balay   PetscScalar       *work,*v,*vv;
601aa05aa95SBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;
6027056b6fcSBarry Smith 
6033a40ed3dSBarry Smith   PetscFunctionBegin;
60439ddd567SLois Curfman McInnes   if (mdn->size == 1) {
60539ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
606aa05aa95SBarry Smith   } else {
607aa05aa95SBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
608ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr);
609ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
610aa05aa95SBarry Smith 
611aa05aa95SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
612f4403165SShri Abhyankar     if (format == PETSC_VIEWER_NATIVE) {
613aa05aa95SBarry Smith 
614aa05aa95SBarry Smith       if (!rank) {
615aa05aa95SBarry Smith         /* store the matrix as a dense matrix */
6160700a824SBarry Smith         header[0] = MAT_FILE_CLASSID;
617d0f46423SBarry Smith         header[1] = mat->rmap->N;
618aa05aa95SBarry Smith         header[2] = N;
619aa05aa95SBarry Smith         header[3] = MATRIX_BINARY_FORMAT_DENSE;
620aa05aa95SBarry Smith         ierr      = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
621aa05aa95SBarry Smith 
622aa05aa95SBarry Smith         /* get largest work array needed for transposing array */
623d0f46423SBarry Smith         mmax = mat->rmap->n;
624aa05aa95SBarry Smith         for (i=1; i<size; i++) {
625d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
6268965ea79SLois Curfman McInnes         }
627785e854fSJed Brown         ierr = PetscMalloc1(mmax*N,&work);CHKERRQ(ierr);
628aa05aa95SBarry Smith 
629aa05aa95SBarry Smith         /* write out local array, by rows */
630d0f46423SBarry Smith         m = mat->rmap->n;
631aa05aa95SBarry Smith         v = a->v;
632aa05aa95SBarry Smith         for (j=0; j<N; j++) {
633aa05aa95SBarry Smith           for (i=0; i<m; i++) {
634578230a0SSatish Balay             work[j + i*N] = *v++;
635aa05aa95SBarry Smith           }
636aa05aa95SBarry Smith         }
637aa05aa95SBarry Smith         ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
638aa05aa95SBarry Smith         /* get largest work array to receive messages from other processes, excludes process zero */
639aa05aa95SBarry Smith         mmax = 0;
640aa05aa95SBarry Smith         for (i=1; i<size; i++) {
641d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
642aa05aa95SBarry Smith         }
643785e854fSJed Brown         ierr = PetscMalloc1(mmax*N,&vv);CHKERRQ(ierr);
644aa05aa95SBarry Smith         for (k = 1; k < size; k++) {
645f8009846SMatthew Knepley           v    = vv;
646d0f46423SBarry Smith           m    = mat->rmap->range[k+1] - mat->rmap->range[k];
647ce94432eSBarry Smith           ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
648aa05aa95SBarry Smith 
649aa05aa95SBarry Smith           for (j = 0; j < N; j++) {
650aa05aa95SBarry Smith             for (i = 0; i < m; i++) {
651578230a0SSatish Balay               work[j + i*N] = *v++;
652aa05aa95SBarry Smith             }
653aa05aa95SBarry Smith           }
654aa05aa95SBarry Smith           ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
655aa05aa95SBarry Smith         }
656aa05aa95SBarry Smith         ierr = PetscFree(work);CHKERRQ(ierr);
657578230a0SSatish Balay         ierr = PetscFree(vv);CHKERRQ(ierr);
658aa05aa95SBarry Smith       } else {
659ce94432eSBarry Smith         ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
660aa05aa95SBarry Smith       }
6616a9046bcSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE)");
662aa05aa95SBarry Smith   }
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
6648965ea79SLois Curfman McInnes }
6658965ea79SLois Curfman McInnes 
6667da1fb6eSBarry Smith extern PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
6679804daf3SBarry Smith #include <petscdraw.h>
6686849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6698965ea79SLois Curfman McInnes {
67039ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
671dfbe8321SBarry Smith   PetscErrorCode    ierr;
6727da1fb6eSBarry Smith   PetscMPIInt       rank = mdn->rank;
67319fd82e9SBarry Smith   PetscViewerType   vtype;
674ace3abfcSBarry Smith   PetscBool         iascii,isdraw;
675b0a32e0cSBarry Smith   PetscViewer       sviewer;
676f3ef73ceSBarry Smith   PetscViewerFormat format;
6778965ea79SLois Curfman McInnes 
6783a40ed3dSBarry Smith   PetscFunctionBegin;
679251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
680251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
68132077d6dSBarry Smith   if (iascii) {
682b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
683b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
684456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6854e220ebcSLois Curfman McInnes       MatInfo info;
686888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
6871575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
6887b23a99aSBarry 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);
689b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6901575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
6915d00a290SHong Zhang       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
6923a40ed3dSBarry Smith       PetscFunctionReturn(0);
693fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6943a40ed3dSBarry Smith       PetscFunctionReturn(0);
6958965ea79SLois Curfman McInnes     }
696f1af5d2fSBarry Smith   } else if (isdraw) {
697b0a32e0cSBarry Smith     PetscDraw draw;
698ace3abfcSBarry Smith     PetscBool isnull;
699f1af5d2fSBarry Smith 
700b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
701b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
702f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
703f1af5d2fSBarry Smith   }
70477ed5343SBarry Smith 
7057da1fb6eSBarry Smith   {
7068965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
7078965ea79SLois Curfman McInnes     Mat         A;
708d0f46423SBarry Smith     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
709ba8c8a56SBarry Smith     PetscInt    *cols;
710ba8c8a56SBarry Smith     PetscScalar *vals;
7118965ea79SLois Curfman McInnes 
712ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr);
7138965ea79SLois Curfman McInnes     if (!rank) {
714f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
7153a40ed3dSBarry Smith     } else {
716f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
7178965ea79SLois Curfman McInnes     }
7187adad957SLisandro Dalcin     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
719878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
7200298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(A,NULL);CHKERRQ(ierr);
7213bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)A);CHKERRQ(ierr);
7228965ea79SLois Curfman McInnes 
72339ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
72439ddd567SLois Curfman McInnes        but it's quick for now */
72551022da4SBarry Smith     A->insertmode = INSERT_VALUES;
7262205254eSKarl Rupp 
7272205254eSKarl Rupp     row = mat->rmap->rstart;
7282205254eSKarl Rupp     m   = mdn->A->rmap->n;
7298965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
730ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
731ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
732ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
73339ddd567SLois Curfman McInnes       row++;
7348965ea79SLois Curfman McInnes     }
7358965ea79SLois Curfman McInnes 
7366d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7376d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7383f08860eSBarry Smith     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
739b9b97703SBarry Smith     if (!rank) {
7401a9d3c3cSBarry Smith       ierr = PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr);
7417da1fb6eSBarry Smith       ierr = MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
7428965ea79SLois Curfman McInnes     }
7433f08860eSBarry Smith     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
744b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7456bf464f9SBarry Smith     ierr = MatDestroy(&A);CHKERRQ(ierr);
7468965ea79SLois Curfman McInnes   }
7473a40ed3dSBarry Smith   PetscFunctionReturn(0);
7488965ea79SLois Curfman McInnes }
7498965ea79SLois Curfman McInnes 
750dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
7518965ea79SLois Curfman McInnes {
752dfbe8321SBarry Smith   PetscErrorCode ierr;
753ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw,issocket;
7548965ea79SLois Curfman McInnes 
755433994e6SBarry Smith   PetscFunctionBegin;
756251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
757251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
758251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr);
759251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
7600f5bd95cSBarry Smith 
76132077d6dSBarry Smith   if (iascii || issocket || isdraw) {
762f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7630f5bd95cSBarry Smith   } else if (isbinary) {
7643a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
76511aeaf0aSBarry Smith   }
7663a40ed3dSBarry Smith   PetscFunctionReturn(0);
7678965ea79SLois Curfman McInnes }
7688965ea79SLois Curfman McInnes 
769dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7708965ea79SLois Curfman McInnes {
7713501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7723501a2bdSLois Curfman McInnes   Mat            mdn  = mat->A;
773dfbe8321SBarry Smith   PetscErrorCode ierr;
774329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7758965ea79SLois Curfman McInnes 
7763a40ed3dSBarry Smith   PetscFunctionBegin;
7774e220ebcSLois Curfman McInnes   info->block_size = 1.0;
7782205254eSKarl Rupp 
7794e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7802205254eSKarl Rupp 
7814e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7824e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7838965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7844e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7854e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7864e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7874e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7884e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7898965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
790b2566f29SBarry Smith     ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
7912205254eSKarl Rupp 
7924e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7934e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7944e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7954e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7964e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7978965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
798b2566f29SBarry Smith     ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
7992205254eSKarl Rupp 
8004e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
8014e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
8024e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
8034e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
8044e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
8058965ea79SLois Curfman McInnes   }
8064e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
8074e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
8084e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
8093a40ed3dSBarry Smith   PetscFunctionReturn(0);
8108965ea79SLois Curfman McInnes }
8118965ea79SLois Curfman McInnes 
812ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg)
8138965ea79SLois Curfman McInnes {
81439ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
815dfbe8321SBarry Smith   PetscErrorCode ierr;
8168965ea79SLois Curfman McInnes 
8173a40ed3dSBarry Smith   PetscFunctionBegin;
81812c028f9SKris Buschelman   switch (op) {
819512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
82012c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
82112c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
82243674050SBarry Smith     MatCheckPreallocated(A,1);
8234e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82412c028f9SKris Buschelman     break;
82512c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
82643674050SBarry Smith     MatCheckPreallocated(A,1);
8274e0d8c25SBarry Smith     a->roworiented = flg;
8284e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82912c028f9SKris Buschelman     break;
8304e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
83113fa8e87SLisandro Dalcin   case MAT_KEEP_NONZERO_PATTERN:
83212c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
833290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
83412c028f9SKris Buschelman     break;
83512c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
8364e0d8c25SBarry Smith     a->donotstash = flg;
83712c028f9SKris Buschelman     break;
83877e54ba9SKris Buschelman   case MAT_SYMMETRIC:
83977e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8409a4540c5SBarry Smith   case MAT_HERMITIAN:
8419a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
842600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
843290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
84477e54ba9SKris Buschelman     break;
84512c028f9SKris Buschelman   default:
846e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
8473a40ed3dSBarry Smith   }
8483a40ed3dSBarry Smith   PetscFunctionReturn(0);
8498965ea79SLois Curfman McInnes }
8508965ea79SLois Curfman McInnes 
8518965ea79SLois Curfman McInnes 
852dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8535b2fa520SLois Curfman McInnes {
8545b2fa520SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)A->data;
8555b2fa520SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)mdn->A->data;
856bca11509SBarry Smith   const PetscScalar *l,*r;
857bca11509SBarry Smith   PetscScalar       x,*v;
858dfbe8321SBarry Smith   PetscErrorCode    ierr;
859d0f46423SBarry Smith   PetscInt          i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n;
8605b2fa520SLois Curfman McInnes 
8615b2fa520SLois Curfman McInnes   PetscFunctionBegin;
86272d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8635b2fa520SLois Curfman McInnes   if (ll) {
86472d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
865e32f2f54SBarry Smith     if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
866bca11509SBarry Smith     ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr);
8675b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8685b2fa520SLois Curfman McInnes       x = l[i];
8695b2fa520SLois Curfman McInnes       v = mat->v + i;
8705b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8715b2fa520SLois Curfman McInnes     }
872bca11509SBarry Smith     ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr);
873efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8745b2fa520SLois Curfman McInnes   }
8755b2fa520SLois Curfman McInnes   if (rr) {
876175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
877e32f2f54SBarry Smith     if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
878ca9f406cSSatish Balay     ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
879ca9f406cSSatish Balay     ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
880bca11509SBarry Smith     ierr = VecGetArrayRead(mdn->lvec,&r);CHKERRQ(ierr);
8815b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8825b2fa520SLois Curfman McInnes       x = r[i];
8835b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8842205254eSKarl Rupp       for (j=0; j<m; j++) (*v++) *= x;
8855b2fa520SLois Curfman McInnes     }
886bca11509SBarry Smith     ierr = VecRestoreArrayRead(mdn->lvec,&r);CHKERRQ(ierr);
887efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8885b2fa520SLois Curfman McInnes   }
8895b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8905b2fa520SLois Curfman McInnes }
8915b2fa520SLois Curfman McInnes 
892dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
893096963f5SLois Curfman McInnes {
8943501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8953501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
896dfbe8321SBarry Smith   PetscErrorCode ierr;
89713f74950SBarry Smith   PetscInt       i,j;
898329f5518SBarry Smith   PetscReal      sum = 0.0;
89987828ca2SBarry Smith   PetscScalar    *v  = mat->v;
9003501a2bdSLois Curfman McInnes 
9013a40ed3dSBarry Smith   PetscFunctionBegin;
9023501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
903064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
9043501a2bdSLois Curfman McInnes   } else {
9053501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
906d0f46423SBarry Smith       for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) {
907329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
9083501a2bdSLois Curfman McInnes       }
909b2566f29SBarry Smith       ierr = MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
9108f1a2a5eSBarry Smith       *nrm = PetscSqrtReal(*nrm);
911dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr);
9123a40ed3dSBarry Smith     } else if (type == NORM_1) {
913329f5518SBarry Smith       PetscReal *tmp,*tmp2;
914dcca6d9dSJed Brown       ierr = PetscMalloc2(A->cmap->N,&tmp,A->cmap->N,&tmp2);CHKERRQ(ierr);
91574ed9c26SBarry Smith       ierr = PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
91674ed9c26SBarry Smith       ierr = PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
917064f8208SBarry Smith       *nrm = 0.0;
9183501a2bdSLois Curfman McInnes       v    = mat->v;
919d0f46423SBarry Smith       for (j=0; j<mdn->A->cmap->n; j++) {
920d0f46423SBarry Smith         for (i=0; i<mdn->A->rmap->n; i++) {
92167e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
9223501a2bdSLois Curfman McInnes         }
9233501a2bdSLois Curfman McInnes       }
924b2566f29SBarry Smith       ierr = MPIU_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
925d0f46423SBarry Smith       for (j=0; j<A->cmap->N; j++) {
926064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
9273501a2bdSLois Curfman McInnes       }
9288627564fSBarry Smith       ierr = PetscFree2(tmp,tmp2);CHKERRQ(ierr);
929d0f46423SBarry Smith       ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
9303a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
931329f5518SBarry Smith       PetscReal ntemp;
9323501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
933b2566f29SBarry Smith       ierr = MPIU_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
934ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm");
9353501a2bdSLois Curfman McInnes   }
9363a40ed3dSBarry Smith   PetscFunctionReturn(0);
9373501a2bdSLois Curfman McInnes }
9383501a2bdSLois Curfman McInnes 
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;
950cf37664fSBarry Smith   if (reuse == MAT_INPLACE_MATRIX  && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place");
951cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
952ce94432eSBarry 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;
961785e854fSJed Brown   ierr = PetscMalloc1(m,&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);
970cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
9713501a2bdSLois Curfman McInnes     *matout = B;
9723501a2bdSLois Curfman McInnes   } else {
97328be2f97SBarry 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 
9824994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A)
983273d9f13SBarry Smith {
984dfbe8321SBarry Smith   PetscErrorCode ierr;
985273d9f13SBarry Smith 
986273d9f13SBarry Smith   PetscFunctionBegin;
987273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
988273d9f13SBarry Smith   PetscFunctionReturn(0);
989273d9f13SBarry Smith }
990273d9f13SBarry Smith 
991488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
992488007eeSBarry Smith {
993488007eeSBarry Smith   PetscErrorCode ierr;
994488007eeSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;
995488007eeSBarry Smith 
996488007eeSBarry Smith   PetscFunctionBegin;
997488007eeSBarry Smith   ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr);
998a3fa217bSJose E. Roman   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
999488007eeSBarry Smith   PetscFunctionReturn(0);
1000488007eeSBarry Smith }
1001488007eeSBarry Smith 
10027087cfbeSBarry Smith PetscErrorCode  MatConjugate_MPIDense(Mat mat)
1003ba337c44SJed Brown {
1004ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)mat->data;
1005ba337c44SJed Brown   PetscErrorCode ierr;
1006ba337c44SJed Brown 
1007ba337c44SJed Brown   PetscFunctionBegin;
1008ba337c44SJed Brown   ierr = MatConjugate(a->A);CHKERRQ(ierr);
1009ba337c44SJed Brown   PetscFunctionReturn(0);
1010ba337c44SJed Brown }
1011ba337c44SJed Brown 
1012ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A)
1013ba337c44SJed Brown {
1014ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1015ba337c44SJed Brown   PetscErrorCode ierr;
1016ba337c44SJed Brown 
1017ba337c44SJed Brown   PetscFunctionBegin;
1018ba337c44SJed Brown   ierr = MatRealPart(a->A);CHKERRQ(ierr);
1019ba337c44SJed Brown   PetscFunctionReturn(0);
1020ba337c44SJed Brown }
1021ba337c44SJed Brown 
1022ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1023ba337c44SJed Brown {
1024ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1025ba337c44SJed Brown   PetscErrorCode ierr;
1026ba337c44SJed Brown 
1027ba337c44SJed Brown   PetscFunctionBegin;
1028ba337c44SJed Brown   ierr = MatImaginaryPart(a->A);CHKERRQ(ierr);
1029ba337c44SJed Brown   PetscFunctionReturn(0);
1030ba337c44SJed Brown }
1031ba337c44SJed Brown 
10320716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);
10330716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
10340716a85fSBarry Smith {
10350716a85fSBarry Smith   PetscErrorCode ierr;
10360716a85fSBarry Smith   PetscInt       i,n;
10370716a85fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
10380716a85fSBarry Smith   PetscReal      *work;
10390716a85fSBarry Smith 
10400716a85fSBarry Smith   PetscFunctionBegin;
10410298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr);
1042785e854fSJed Brown   ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
10430716a85fSBarry Smith   ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr);
10440716a85fSBarry Smith   if (type == NORM_2) {
10450716a85fSBarry Smith     for (i=0; i<n; i++) work[i] *= work[i];
10460716a85fSBarry Smith   }
10470716a85fSBarry Smith   if (type == NORM_INFINITY) {
1048b2566f29SBarry Smith     ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr);
10490716a85fSBarry Smith   } else {
1050b2566f29SBarry Smith     ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr);
10510716a85fSBarry Smith   }
10520716a85fSBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
10530716a85fSBarry Smith   if (type == NORM_2) {
10548f1a2a5eSBarry Smith     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
10550716a85fSBarry Smith   }
10560716a85fSBarry Smith   PetscFunctionReturn(0);
10570716a85fSBarry Smith }
10580716a85fSBarry Smith 
105973a71a0fSBarry Smith static PetscErrorCode  MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
106073a71a0fSBarry Smith {
106173a71a0fSBarry Smith   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;
106273a71a0fSBarry Smith   PetscErrorCode ierr;
106373a71a0fSBarry Smith   PetscScalar    *a;
106473a71a0fSBarry Smith   PetscInt       m,n,i;
106573a71a0fSBarry Smith 
106673a71a0fSBarry Smith   PetscFunctionBegin;
106773a71a0fSBarry Smith   ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr);
10688c778c55SBarry Smith   ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr);
106973a71a0fSBarry Smith   for (i=0; i<m*n; i++) {
107073a71a0fSBarry Smith     ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr);
107173a71a0fSBarry Smith   }
10728c778c55SBarry Smith   ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr);
107373a71a0fSBarry Smith   PetscFunctionReturn(0);
107473a71a0fSBarry Smith }
107573a71a0fSBarry Smith 
1076fd4e9aacSBarry Smith extern PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat,Mat);
1077fd4e9aacSBarry Smith 
10783b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MPIDense(Mat A,PetscBool  *missing,PetscInt *d)
10793b49f96aSBarry Smith {
10803b49f96aSBarry Smith   PetscFunctionBegin;
10813b49f96aSBarry Smith   *missing = PETSC_FALSE;
10823b49f96aSBarry Smith   PetscFunctionReturn(0);
10833b49f96aSBarry Smith }
10843b49f96aSBarry Smith 
10858965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
108609dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense,
108709dc0095SBarry Smith                                         MatGetRow_MPIDense,
108809dc0095SBarry Smith                                         MatRestoreRow_MPIDense,
108909dc0095SBarry Smith                                         MatMult_MPIDense,
109097304618SKris Buschelman                                 /*  4*/ MatMultAdd_MPIDense,
10917c922b88SBarry Smith                                         MatMultTranspose_MPIDense,
10927c922b88SBarry Smith                                         MatMultTransposeAdd_MPIDense,
10938965ea79SLois Curfman McInnes                                         0,
109409dc0095SBarry Smith                                         0,
109509dc0095SBarry Smith                                         0,
109697304618SKris Buschelman                                 /* 10*/ 0,
109709dc0095SBarry Smith                                         0,
109809dc0095SBarry Smith                                         0,
109909dc0095SBarry Smith                                         0,
110009dc0095SBarry Smith                                         MatTranspose_MPIDense,
110197304618SKris Buschelman                                 /* 15*/ MatGetInfo_MPIDense,
11026e4ee0c6SHong Zhang                                         MatEqual_MPIDense,
110309dc0095SBarry Smith                                         MatGetDiagonal_MPIDense,
11045b2fa520SLois Curfman McInnes                                         MatDiagonalScale_MPIDense,
110509dc0095SBarry Smith                                         MatNorm_MPIDense,
110697304618SKris Buschelman                                 /* 20*/ MatAssemblyBegin_MPIDense,
110709dc0095SBarry Smith                                         MatAssemblyEnd_MPIDense,
110809dc0095SBarry Smith                                         MatSetOption_MPIDense,
110909dc0095SBarry Smith                                         MatZeroEntries_MPIDense,
1110d519adbfSMatthew Knepley                                 /* 24*/ MatZeroRows_MPIDense,
1111919b68f7SBarry Smith                                         0,
111201b82886SBarry Smith                                         0,
111301b82886SBarry Smith                                         0,
111401b82886SBarry Smith                                         0,
11154994cf47SJed Brown                                 /* 29*/ MatSetUp_MPIDense,
1116273d9f13SBarry Smith                                         0,
111709dc0095SBarry Smith                                         0,
1118c56a70eeSBarry Smith                                         MatGetDiagonalBlock_MPIDense,
11198c778c55SBarry Smith                                         0,
1120d519adbfSMatthew Knepley                                 /* 34*/ MatDuplicate_MPIDense,
112109dc0095SBarry Smith                                         0,
112209dc0095SBarry Smith                                         0,
112309dc0095SBarry Smith                                         0,
112409dc0095SBarry Smith                                         0,
1125d519adbfSMatthew Knepley                                 /* 39*/ MatAXPY_MPIDense,
11267dae84e0SHong Zhang                                         MatCreateSubMatrices_MPIDense,
112709dc0095SBarry Smith                                         0,
112809dc0095SBarry Smith                                         MatGetValues_MPIDense,
112909dc0095SBarry Smith                                         0,
1130d519adbfSMatthew Knepley                                 /* 44*/ 0,
113109dc0095SBarry Smith                                         MatScale_MPIDense,
11327d68702bSBarry Smith                                         MatShift_Basic,
113309dc0095SBarry Smith                                         0,
113409dc0095SBarry Smith                                         0,
113573a71a0fSBarry Smith                                 /* 49*/ MatSetRandom_MPIDense,
113609dc0095SBarry Smith                                         0,
113709dc0095SBarry Smith                                         0,
113809dc0095SBarry Smith                                         0,
113909dc0095SBarry Smith                                         0,
1140d519adbfSMatthew Knepley                                 /* 54*/ 0,
114109dc0095SBarry Smith                                         0,
114209dc0095SBarry Smith                                         0,
114309dc0095SBarry Smith                                         0,
114409dc0095SBarry Smith                                         0,
11457dae84e0SHong Zhang                                 /* 59*/ MatCreateSubMatrix_MPIDense,
1146b9b97703SBarry Smith                                         MatDestroy_MPIDense,
1147b9b97703SBarry Smith                                         MatView_MPIDense,
1148357abbc8SBarry Smith                                         0,
114997304618SKris Buschelman                                         0,
1150d519adbfSMatthew Knepley                                 /* 64*/ 0,
115197304618SKris Buschelman                                         0,
115297304618SKris Buschelman                                         0,
115397304618SKris Buschelman                                         0,
115497304618SKris Buschelman                                         0,
1155d519adbfSMatthew Knepley                                 /* 69*/ 0,
115697304618SKris Buschelman                                         0,
115797304618SKris Buschelman                                         0,
115897304618SKris Buschelman                                         0,
115997304618SKris Buschelman                                         0,
1160d519adbfSMatthew Knepley                                 /* 74*/ 0,
116197304618SKris Buschelman                                         0,
116297304618SKris Buschelman                                         0,
116397304618SKris Buschelman                                         0,
116497304618SKris Buschelman                                         0,
1165d519adbfSMatthew Knepley                                 /* 79*/ 0,
116697304618SKris Buschelman                                         0,
116797304618SKris Buschelman                                         0,
116897304618SKris Buschelman                                         0,
11695bba2384SShri Abhyankar                                 /* 83*/ MatLoad_MPIDense,
1170865e5f61SKris Buschelman                                         0,
1171865e5f61SKris Buschelman                                         0,
1172865e5f61SKris Buschelman                                         0,
1173865e5f61SKris Buschelman                                         0,
1174865e5f61SKris Buschelman                                         0,
11759775878aSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
1176320f2790SHong Zhang                                 /* 89*/ MatMatMult_MPIDense_MPIDense,
1177320f2790SHong Zhang                                         MatMatMultSymbolic_MPIDense_MPIDense,
11789775878aSHong Zhang #else
11799775878aSHong Zhang                                 /* 89*/ 0,
1180865e5f61SKris Buschelman                                         0,
11819775878aSHong Zhang #endif
1182fd4e9aacSBarry Smith                                         MatMatMultNumeric_MPIDense,
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,
12043b49f96aSBarry Smith                                         MatMissingDiagonal_MPIDense,
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,
12163964eb88SJed Brown                                         MatGetColumnNorms_MPIDense,
12173964eb88SJed Brown                                         0,
12183964eb88SJed Brown                                         0,
12193964eb88SJed Brown                                         0,
12203964eb88SJed Brown                                 /*129*/ 0,
1221cb20be35SHong Zhang                                         MatTransposeMatMult_MPIDense_MPIDense,
1222cb20be35SHong Zhang                                         MatTransposeMatMultSymbolic_MPIDense_MPIDense,
1223cb20be35SHong Zhang                                         MatTransposeMatMultNumeric_MPIDense_MPIDense,
12243964eb88SJed Brown                                         0,
12253964eb88SJed Brown                                 /*134*/ 0,
12263964eb88SJed Brown                                         0,
12273964eb88SJed Brown                                         0,
12283964eb88SJed Brown                                         0,
12293964eb88SJed Brown                                         0,
12303964eb88SJed Brown                                 /*139*/ 0,
1231f9426fe0SMark Adams                                         0,
12323964eb88SJed Brown                                         0
1233ba337c44SJed Brown };
12348965ea79SLois Curfman McInnes 
12357087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1236a23d5eceSKris Buschelman {
1237a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1238dfbe8321SBarry Smith   PetscErrorCode ierr;
1239a23d5eceSKris Buschelman 
1240a23d5eceSKris Buschelman   PetscFunctionBegin;
1241a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1242a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1243a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1244a23d5eceSKris Buschelman 
1245a23d5eceSKris Buschelman   a       = (Mat_MPIDense*)mat->data;
124634ef9618SShri Abhyankar   ierr    = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
124734ef9618SShri Abhyankar   ierr    = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
124834ef9618SShri Abhyankar   a->nvec = mat->cmap->n;
124934ef9618SShri Abhyankar 
1250f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1251d0f46423SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr);
12525c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
12535c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
12543bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr);
1255a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1256a23d5eceSKris Buschelman }
1257a23d5eceSKris Buschelman 
125865b80a83SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
1259cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
12608baccfbdSHong Zhang {
12618ea901baSHong Zhang   Mat            mat_elemental;
12628ea901baSHong Zhang   PetscErrorCode ierr;
126332d7a744SHong Zhang   PetscScalar    *v;
126432d7a744SHong Zhang   PetscInt       m=A->rmap->n,N=A->cmap->N,rstart=A->rmap->rstart,i,*rows,*cols;
12658ea901baSHong Zhang 
12668baccfbdSHong Zhang   PetscFunctionBegin;
1267378336b6SHong Zhang   if (reuse == MAT_REUSE_MATRIX) {
1268378336b6SHong Zhang     mat_elemental = *newmat;
1269378336b6SHong Zhang     ierr = MatZeroEntries(*newmat);CHKERRQ(ierr);
1270378336b6SHong Zhang   } else {
1271378336b6SHong Zhang     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
1272378336b6SHong Zhang     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1273378336b6SHong Zhang     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
1274378336b6SHong Zhang     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
127532d7a744SHong Zhang     ierr = MatSetOption(mat_elemental,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
1276378336b6SHong Zhang   }
1277378336b6SHong Zhang 
127832d7a744SHong Zhang   ierr = PetscMalloc2(m,&rows,N,&cols);CHKERRQ(ierr);
127932d7a744SHong Zhang   for (i=0; i<N; i++) cols[i] = i;
128032d7a744SHong Zhang   for (i=0; i<m; i++) rows[i] = rstart + i;
12818ea901baSHong Zhang 
12828ea901baSHong Zhang   /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
128332d7a744SHong Zhang   ierr = MatDenseGetArray(A,&v);CHKERRQ(ierr);
128432d7a744SHong Zhang   ierr = MatSetValues(mat_elemental,m,rows,N,cols,v,ADD_VALUES);CHKERRQ(ierr);
12858ea901baSHong Zhang   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12868ea901baSHong Zhang   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
128732d7a744SHong Zhang   ierr = MatDenseRestoreArray(A,&v);CHKERRQ(ierr);
128832d7a744SHong Zhang   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
12898ea901baSHong Zhang 
1290511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
129128be2f97SBarry Smith     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
12928ea901baSHong Zhang   } else {
12938ea901baSHong Zhang     *newmat = mat_elemental;
12948ea901baSHong Zhang   }
12958baccfbdSHong Zhang   PetscFunctionReturn(0);
12968baccfbdSHong Zhang }
129765b80a83SHong Zhang #endif
12988baccfbdSHong Zhang 
12998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat)
1300273d9f13SBarry Smith {
1301273d9f13SBarry Smith   Mat_MPIDense   *a;
1302dfbe8321SBarry Smith   PetscErrorCode ierr;
1303273d9f13SBarry Smith 
1304273d9f13SBarry Smith   PetscFunctionBegin;
1305b00a9115SJed Brown   ierr      = PetscNewLog(mat,&a);CHKERRQ(ierr);
1306b0a32e0cSBarry Smith   mat->data = (void*)a;
1307273d9f13SBarry Smith   ierr      = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1308273d9f13SBarry Smith 
1309273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1310ce94432eSBarry Smith   ierr            = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);CHKERRQ(ierr);
1311ce94432eSBarry Smith   ierr            = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);CHKERRQ(ierr);
1312273d9f13SBarry Smith 
1313273d9f13SBarry Smith   /* build cache for off array entries formed */
1314273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
13152205254eSKarl Rupp 
1316ce94432eSBarry Smith   ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);CHKERRQ(ierr);
1317273d9f13SBarry Smith 
1318273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1319273d9f13SBarry Smith   a->lvec        = 0;
1320273d9f13SBarry Smith   a->Mvctx       = 0;
1321273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1322273d9f13SBarry Smith 
1323bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);CHKERRQ(ierr);
1324*8572280aSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr);
1325*8572280aSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",MatDenseGetArrayRead_MPIDense);CHKERRQ(ierr);
1326*8572280aSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",MatDenseRestoreArrayRead_MPIDense);CHKERRQ(ierr);
1327d3042a70SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",MatDensePlaceArray_MPIDense);CHKERRQ(ierr);
1328d3042a70SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",MatDenseResetArray_MPIDense);CHKERRQ(ierr);
13298baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
13308baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",MatConvert_MPIDense_Elemental);CHKERRQ(ierr);
13318baccfbdSHong Zhang #endif
1332bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1333bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
1334bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
1335bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
13368949adfdSHong Zhang 
13378949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",MatTransposeMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
13388949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",MatTransposeMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
13398949adfdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",MatTransposeMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
134038aed534SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr);
1341273d9f13SBarry Smith   PetscFunctionReturn(0);
1342273d9f13SBarry Smith }
1343273d9f13SBarry Smith 
1344209238afSKris Buschelman /*MC
1345002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1346209238afSKris Buschelman 
1347209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1348209238afSKris Buschelman    and MATMPIDENSE otherwise.
1349209238afSKris Buschelman 
1350209238afSKris Buschelman    Options Database Keys:
1351209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1352209238afSKris Buschelman 
1353209238afSKris Buschelman   Level: beginner
1354209238afSKris Buschelman 
135501b82886SBarry Smith 
1356209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1357209238afSKris Buschelman M*/
1358209238afSKris Buschelman 
1359273d9f13SBarry Smith /*@C
1360273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1361273d9f13SBarry Smith 
1362273d9f13SBarry Smith    Not collective
1363273d9f13SBarry Smith 
1364273d9f13SBarry Smith    Input Parameters:
13651c4f3114SJed Brown .  B - the matrix
13660298fd71SBarry Smith -  data - optional location of matrix data.  Set data=NULL for PETSc
1367273d9f13SBarry Smith    to control all matrix memory allocation.
1368273d9f13SBarry Smith 
1369273d9f13SBarry Smith    Notes:
1370273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1371273d9f13SBarry Smith    storage by columns.
1372273d9f13SBarry Smith 
1373273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1374273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
13750298fd71SBarry Smith    set data=NULL.
1376273d9f13SBarry Smith 
1377273d9f13SBarry Smith    Level: intermediate
1378273d9f13SBarry Smith 
1379273d9f13SBarry Smith .keywords: matrix,dense, parallel
1380273d9f13SBarry Smith 
1381273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1382273d9f13SBarry Smith @*/
13831c4f3114SJed Brown PetscErrorCode  MatMPIDenseSetPreallocation(Mat B,PetscScalar *data)
1384273d9f13SBarry Smith {
13854ac538c5SBarry Smith   PetscErrorCode ierr;
1386273d9f13SBarry Smith 
1387273d9f13SBarry Smith   PetscFunctionBegin;
13881c4f3114SJed Brown   ierr = PetscTryMethod(B,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(B,data));CHKERRQ(ierr);
1389273d9f13SBarry Smith   PetscFunctionReturn(0);
1390273d9f13SBarry Smith }
1391273d9f13SBarry Smith 
1392d3042a70SBarry Smith /*@
1393d3042a70SBarry Smith    MatDensePlaceArray - Allows one to replace the array in a dense array with an
1394d3042a70SBarry Smith    array provided by the user. This is useful to avoid copying an array
1395d3042a70SBarry Smith    into a matrix
1396d3042a70SBarry Smith 
1397d3042a70SBarry Smith    Not Collective
1398d3042a70SBarry Smith 
1399d3042a70SBarry Smith    Input Parameters:
1400d3042a70SBarry Smith +  mat - the matrix
1401d3042a70SBarry Smith -  array - the array in column major order
1402d3042a70SBarry Smith 
1403d3042a70SBarry Smith    Notes:
1404d3042a70SBarry Smith    You can return to the original array with a call to MatDenseResetArray(). The user is responsible for freeing this array; it will not be
1405d3042a70SBarry Smith    freed when the matrix is destroyed.
1406d3042a70SBarry Smith 
1407d3042a70SBarry Smith    Level: developer
1408d3042a70SBarry Smith 
1409d3042a70SBarry Smith .seealso: MatDenseGetArray(), MatDenseResetArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()
1410d3042a70SBarry Smith 
1411d3042a70SBarry Smith @*/
1412d3042a70SBarry Smith PetscErrorCode  MatDensePlaceArray(Mat mat,const PetscScalar array[])
1413d3042a70SBarry Smith {
1414d3042a70SBarry Smith   PetscErrorCode ierr;
1415d3042a70SBarry Smith   PetscFunctionBegin;
1416d3042a70SBarry Smith   ierr = PetscUseMethod(mat,"MatDensePlaceArray_C",(Mat,const PetscScalar*),(mat,array));CHKERRQ(ierr);
1417d3042a70SBarry Smith   ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr);
1418d3042a70SBarry Smith   PetscFunctionReturn(0);
1419d3042a70SBarry Smith }
1420d3042a70SBarry Smith 
1421d3042a70SBarry Smith /*@
1422d3042a70SBarry Smith    MatDenseResetArray - Resets the matrix array to that it previously had before the call to MatDensePlaceArray()
1423d3042a70SBarry Smith 
1424d3042a70SBarry Smith    Not Collective
1425d3042a70SBarry Smith 
1426d3042a70SBarry Smith    Input Parameters:
1427d3042a70SBarry Smith .  mat - the matrix
1428d3042a70SBarry Smith 
1429d3042a70SBarry Smith    Notes:
1430d3042a70SBarry Smith    You can only call this after a call to MatDensePlaceArray()
1431d3042a70SBarry Smith 
1432d3042a70SBarry Smith    Level: developer
1433d3042a70SBarry Smith 
1434d3042a70SBarry Smith .seealso: MatDenseGetArray(), MatDensePlaceArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()
1435d3042a70SBarry Smith 
1436d3042a70SBarry Smith @*/
1437d3042a70SBarry Smith PetscErrorCode  MatDenseResetArray(Mat mat)
1438d3042a70SBarry Smith {
1439d3042a70SBarry Smith   PetscErrorCode ierr;
1440d3042a70SBarry Smith   PetscFunctionBegin;
1441d3042a70SBarry Smith   ierr = PetscUseMethod(mat,"MatDenseResetArray_C",(Mat),(mat));CHKERRQ(ierr);
1442d3042a70SBarry Smith   ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr);
1443d3042a70SBarry Smith   PetscFunctionReturn(0);
1444d3042a70SBarry Smith }
1445d3042a70SBarry Smith 
14468965ea79SLois Curfman McInnes /*@C
144769b1f4b7SBarry Smith    MatCreateDense - Creates a parallel matrix in dense format.
14488965ea79SLois Curfman McInnes 
1449db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1450db81eaa0SLois Curfman McInnes 
14518965ea79SLois Curfman McInnes    Input Parameters:
1452db81eaa0SLois Curfman McInnes +  comm - MPI communicator
14538965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1454db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
14558965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1456db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
14576cfe35ddSJose E. Roman -  data - optional location of matrix data.  Set data=NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1458dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
14598965ea79SLois Curfman McInnes 
14608965ea79SLois Curfman McInnes    Output Parameter:
1461477f1c0bSLois Curfman McInnes .  A - the matrix
14628965ea79SLois Curfman McInnes 
1463b259b22eSLois Curfman McInnes    Notes:
146439ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
146539ddd567SLois Curfman McInnes    storage by columns.
14668965ea79SLois Curfman McInnes 
146718f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
146818f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
14696cfe35ddSJose E. Roman    set data=NULL (PETSC_NULL_SCALAR for Fortran users).
147018f449edSLois Curfman McInnes 
14718965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
14728965ea79SLois Curfman McInnes    (possibly both).
14738965ea79SLois Curfman McInnes 
1474027ccd11SLois Curfman McInnes    Level: intermediate
1475027ccd11SLois Curfman McInnes 
147639ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
14778965ea79SLois Curfman McInnes 
147839ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
14798965ea79SLois Curfman McInnes @*/
148069b1f4b7SBarry Smith PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
14818965ea79SLois Curfman McInnes {
14826849ba73SBarry Smith   PetscErrorCode ierr;
148313f74950SBarry Smith   PetscMPIInt    size;
14848965ea79SLois Curfman McInnes 
14853a40ed3dSBarry Smith   PetscFunctionBegin;
1486f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1487f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1488273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1489273d9f13SBarry Smith   if (size > 1) {
1490273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1491273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
14926cfe35ddSJose E. Roman     if (data) {  /* user provided data array, so no need to assemble */
14936cfe35ddSJose E. Roman       ierr = MatSetUpMultiply_MPIDense(*A);CHKERRQ(ierr);
14946cfe35ddSJose E. Roman       (*A)->assembled = PETSC_TRUE;
14956cfe35ddSJose E. Roman     }
1496273d9f13SBarry Smith   } else {
1497273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1498273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
14998c469469SLois Curfman McInnes   }
15003a40ed3dSBarry Smith   PetscFunctionReturn(0);
15018965ea79SLois Curfman McInnes }
15028965ea79SLois Curfman McInnes 
15036849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
15048965ea79SLois Curfman McInnes {
15058965ea79SLois Curfman McInnes   Mat            mat;
15063501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1507dfbe8321SBarry Smith   PetscErrorCode ierr;
15088965ea79SLois Curfman McInnes 
15093a40ed3dSBarry Smith   PetscFunctionBegin;
15108965ea79SLois Curfman McInnes   *newmat = 0;
1511ce94432eSBarry Smith   ierr    = MatCreate(PetscObjectComm((PetscObject)A),&mat);CHKERRQ(ierr);
1512d0f46423SBarry Smith   ierr    = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
15137adad957SLisandro Dalcin   ierr    = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1514834f8fabSBarry Smith   a       = (Mat_MPIDense*)mat->data;
1515e04c1aa4SHong Zhang   ierr    = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
15165aa7edbeSHong Zhang 
1517d5f3da31SBarry Smith   mat->factortype   = A->factortype;
1518c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1519273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
15208965ea79SLois Curfman McInnes 
15218965ea79SLois Curfman McInnes   a->size         = oldmat->size;
15228965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1523e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1524b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
15253782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1526e04c1aa4SHong Zhang 
15271e1e43feSBarry Smith   ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr);
15281e1e43feSBarry Smith   ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr);
15298965ea79SLois Curfman McInnes 
1530329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
15315609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
15323bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr);
153301b82886SBarry Smith 
15348965ea79SLois Curfman McInnes   *newmat = mat;
15353a40ed3dSBarry Smith   PetscFunctionReturn(0);
15368965ea79SLois Curfman McInnes }
15378965ea79SLois Curfman McInnes 
15389d36ed5fSBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat)
153986d161a7SShri Abhyankar {
154086d161a7SShri Abhyankar   PetscErrorCode ierr;
154186d161a7SShri Abhyankar   PetscMPIInt    rank,size;
154274c13d95SBarry Smith   const PetscInt *rowners;
154374c13d95SBarry Smith   PetscInt       i,m,n,nz,j,mMax;
154486d161a7SShri Abhyankar   PetscScalar    *array,*vals,*vals_ptr;
15459d36ed5fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)newmat->data;
154686d161a7SShri Abhyankar 
154786d161a7SShri Abhyankar   PetscFunctionBegin;
154886d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
154986d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
155086d161a7SShri Abhyankar 
155174c13d95SBarry Smith   /* determine ownership of rows and columns */
155274c13d95SBarry Smith   m = (newmat->rmap->n < 0) ? PETSC_DECIDE : newmat->rmap->n;
155374c13d95SBarry Smith   n = (newmat->cmap->n < 0) ? PETSC_DECIDE : newmat->cmap->n;
155486d161a7SShri Abhyankar 
155574c13d95SBarry Smith   ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr);
15569d36ed5fSBarry Smith   if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) {
15570298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
15589d36ed5fSBarry Smith   }
15598c778c55SBarry Smith   ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr);
156049467e7dSBarry Smith   ierr = MatGetLocalSize(newmat,&m,NULL);CHKERRQ(ierr);
156174c13d95SBarry Smith   ierr = MatGetOwnershipRanges(newmat,&rowners);CHKERRQ(ierr);
1562396e826eSBarry Smith   ierr = MPI_Reduce(&m,&mMax,1,MPIU_INT,MPI_MAX,0,comm);CHKERRQ(ierr);
156386d161a7SShri Abhyankar   if (!rank) {
15649d36ed5fSBarry Smith     ierr = PetscMalloc1(mMax*N,&vals);CHKERRQ(ierr);
156586d161a7SShri Abhyankar 
156686d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
156786d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
156886d161a7SShri Abhyankar 
156986d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
157086d161a7SShri Abhyankar     vals_ptr = vals;
157186d161a7SShri Abhyankar     for (i=0; i<m; i++) {
157286d161a7SShri Abhyankar       for (j=0; j<N; j++) {
157386d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
157486d161a7SShri Abhyankar       }
157586d161a7SShri Abhyankar     }
157686d161a7SShri Abhyankar 
157786d161a7SShri Abhyankar     /* read in other processors and ship out */
157886d161a7SShri Abhyankar     for (i=1; i<size; i++) {
157986d161a7SShri Abhyankar       nz   = (rowners[i+1] - rowners[i])*N;
158086d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1581a25532f0SBarry Smith       ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
158286d161a7SShri Abhyankar     }
158386d161a7SShri Abhyankar   } else {
158486d161a7SShri Abhyankar     /* receive numeric values */
1585785e854fSJed Brown     ierr = PetscMalloc1(m*N,&vals);CHKERRQ(ierr);
158686d161a7SShri Abhyankar 
158786d161a7SShri Abhyankar     /* receive message of values*/
1588a25532f0SBarry Smith     ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
158986d161a7SShri Abhyankar 
159086d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
159186d161a7SShri Abhyankar     vals_ptr = vals;
159286d161a7SShri Abhyankar     for (i=0; i<m; i++) {
159386d161a7SShri Abhyankar       for (j=0; j<N; j++) {
159486d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
159586d161a7SShri Abhyankar       }
159686d161a7SShri Abhyankar     }
159786d161a7SShri Abhyankar   }
15988c778c55SBarry Smith   ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr);
159986d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
160086d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160186d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160286d161a7SShri Abhyankar   PetscFunctionReturn(0);
160386d161a7SShri Abhyankar }
160486d161a7SShri Abhyankar 
1605112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer)
160686d161a7SShri Abhyankar {
1607dfdea288SBarry Smith   Mat_MPIDense   *a;
160886d161a7SShri Abhyankar   PetscScalar    *vals,*svals;
1609ce94432eSBarry Smith   MPI_Comm       comm;
161086d161a7SShri Abhyankar   MPI_Status     status;
161149467e7dSBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,n,maxnz;
161286d161a7SShri Abhyankar   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
161349467e7dSBarry Smith   PetscInt       *ourlens,*procsnz = 0,jj,*mycols,*smycols;
16149d36ed5fSBarry Smith   PetscInt       i,nz,j,rstart,rend;
161586d161a7SShri Abhyankar   int            fd;
161686d161a7SShri Abhyankar   PetscErrorCode ierr;
161786d161a7SShri Abhyankar 
161886d161a7SShri Abhyankar   PetscFunctionBegin;
1619c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
1620c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
1621ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
162286d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
162386d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
162486d161a7SShri Abhyankar   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
16255872f025SBarry Smith   if (!rank) {
162686d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,(char*)header,4,PETSC_INT);CHKERRQ(ierr);
162786d161a7SShri Abhyankar     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
162886d161a7SShri Abhyankar   }
162986d161a7SShri Abhyankar   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
163086d161a7SShri Abhyankar   M    = header[1]; N = header[2]; nz = header[3];
163186d161a7SShri Abhyankar 
163286d161a7SShri Abhyankar   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
16339d36ed5fSBarry Smith   if (newmat->rmap->N < 0) newmat->rmap->N = M;
16349d36ed5fSBarry Smith   if (newmat->cmap->N < 0) newmat->cmap->N = N;
163586d161a7SShri Abhyankar 
16369d36ed5fSBarry Smith   if (newmat->rmap->N != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%D) and input matrix has (%D)",M,newmat->rmap->N);
16379d36ed5fSBarry Smith   if (newmat->cmap->N != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%D) and input matrix has (%D)",N,newmat->cmap->N);
163886d161a7SShri Abhyankar 
163986d161a7SShri Abhyankar   /*
164086d161a7SShri Abhyankar        Handle case where matrix is stored on disk as a dense matrix
164186d161a7SShri Abhyankar   */
164286d161a7SShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
16439d36ed5fSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
164486d161a7SShri Abhyankar     PetscFunctionReturn(0);
164586d161a7SShri Abhyankar   }
164686d161a7SShri Abhyankar 
164786d161a7SShri Abhyankar   /* determine ownership of all rows */
16482205254eSKarl Rupp   if (newmat->rmap->n < 0) {
16492205254eSKarl Rupp     ierr = PetscMPIIntCast(M/size + ((M % size) > rank),&m);CHKERRQ(ierr);
16502205254eSKarl Rupp   } else {
16512205254eSKarl Rupp     ierr = PetscMPIIntCast(newmat->rmap->n,&m);CHKERRQ(ierr);
16522205254eSKarl Rupp   }
165349467e7dSBarry Smith   if (newmat->cmap->n < 0) {
165449467e7dSBarry Smith     n = PETSC_DECIDE;
165549467e7dSBarry Smith   } else {
165649467e7dSBarry Smith     ierr = PetscMPIIntCast(newmat->cmap->n,&n);CHKERRQ(ierr);
165749467e7dSBarry Smith   }
165849467e7dSBarry Smith 
1659854ce69bSBarry Smith   ierr       = PetscMalloc1(size+2,&rowners);CHKERRQ(ierr);
166086d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
166186d161a7SShri Abhyankar   rowners[0] = 0;
166286d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
166386d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
166486d161a7SShri Abhyankar   }
166586d161a7SShri Abhyankar   rstart = rowners[rank];
166686d161a7SShri Abhyankar   rend   = rowners[rank+1];
166786d161a7SShri Abhyankar 
166886d161a7SShri Abhyankar   /* distribute row lengths to all processors */
166949467e7dSBarry Smith   ierr = PetscMalloc1(rend-rstart,&ourlens);CHKERRQ(ierr);
167086d161a7SShri Abhyankar   if (!rank) {
1671785e854fSJed Brown     ierr = PetscMalloc1(M,&rowlengths);CHKERRQ(ierr);
167286d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1673785e854fSJed Brown     ierr = PetscMalloc1(size,&sndcounts);CHKERRQ(ierr);
167486d161a7SShri Abhyankar     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
167586d161a7SShri Abhyankar     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
167686d161a7SShri Abhyankar     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
167786d161a7SShri Abhyankar   } else {
167886d161a7SShri Abhyankar     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
167986d161a7SShri Abhyankar   }
168086d161a7SShri Abhyankar 
168186d161a7SShri Abhyankar   if (!rank) {
168286d161a7SShri Abhyankar     /* calculate the number of nonzeros on each processor */
1683785e854fSJed Brown     ierr = PetscMalloc1(size,&procsnz);CHKERRQ(ierr);
168486d161a7SShri Abhyankar     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
168586d161a7SShri Abhyankar     for (i=0; i<size; i++) {
168686d161a7SShri Abhyankar       for (j=rowners[i]; j< rowners[i+1]; j++) {
168786d161a7SShri Abhyankar         procsnz[i] += rowlengths[j];
168886d161a7SShri Abhyankar       }
168986d161a7SShri Abhyankar     }
169086d161a7SShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
169186d161a7SShri Abhyankar 
169286d161a7SShri Abhyankar     /* determine max buffer needed and allocate it */
169386d161a7SShri Abhyankar     maxnz = 0;
169486d161a7SShri Abhyankar     for (i=0; i<size; i++) {
169586d161a7SShri Abhyankar       maxnz = PetscMax(maxnz,procsnz[i]);
169686d161a7SShri Abhyankar     }
1697785e854fSJed Brown     ierr = PetscMalloc1(maxnz,&cols);CHKERRQ(ierr);
169886d161a7SShri Abhyankar 
169986d161a7SShri Abhyankar     /* read in my part of the matrix column indices  */
170086d161a7SShri Abhyankar     nz   = procsnz[0];
1701785e854fSJed Brown     ierr = PetscMalloc1(nz,&mycols);CHKERRQ(ierr);
170286d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
170386d161a7SShri Abhyankar 
170486d161a7SShri Abhyankar     /* read in every one elses and ship off */
170586d161a7SShri Abhyankar     for (i=1; i<size; i++) {
170686d161a7SShri Abhyankar       nz   = procsnz[i];
170786d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
170886d161a7SShri Abhyankar       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
170986d161a7SShri Abhyankar     }
171086d161a7SShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
171186d161a7SShri Abhyankar   } else {
171286d161a7SShri Abhyankar     /* determine buffer space needed for message */
171386d161a7SShri Abhyankar     nz = 0;
171486d161a7SShri Abhyankar     for (i=0; i<m; i++) {
171586d161a7SShri Abhyankar       nz += ourlens[i];
171686d161a7SShri Abhyankar     }
1717854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&mycols);CHKERRQ(ierr);
171886d161a7SShri Abhyankar 
171986d161a7SShri Abhyankar     /* receive message of column indices*/
172086d161a7SShri Abhyankar     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
172186d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
172286d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
172386d161a7SShri Abhyankar   }
172486d161a7SShri Abhyankar 
172549467e7dSBarry Smith   ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr);
1726dfdea288SBarry Smith   a = (Mat_MPIDense*)newmat->data;
17272e3566b0SBarry Smith   if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) {
17280298fd71SBarry Smith     ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
1729dfdea288SBarry Smith   }
173086d161a7SShri Abhyankar 
173186d161a7SShri Abhyankar   if (!rank) {
1732785e854fSJed Brown     ierr = PetscMalloc1(maxnz,&vals);CHKERRQ(ierr);
173386d161a7SShri Abhyankar 
173486d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
173586d161a7SShri Abhyankar     nz   = procsnz[0];
173686d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
173786d161a7SShri Abhyankar 
173886d161a7SShri Abhyankar     /* insert into matrix */
173986d161a7SShri Abhyankar     jj      = rstart;
174086d161a7SShri Abhyankar     smycols = mycols;
174186d161a7SShri Abhyankar     svals   = vals;
174286d161a7SShri Abhyankar     for (i=0; i<m; i++) {
174386d161a7SShri Abhyankar       ierr     = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
174486d161a7SShri Abhyankar       smycols += ourlens[i];
174586d161a7SShri Abhyankar       svals   += ourlens[i];
174686d161a7SShri Abhyankar       jj++;
174786d161a7SShri Abhyankar     }
174886d161a7SShri Abhyankar 
174986d161a7SShri Abhyankar     /* read in other processors and ship out */
175086d161a7SShri Abhyankar     for (i=1; i<size; i++) {
175186d161a7SShri Abhyankar       nz   = procsnz[i];
175286d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
175386d161a7SShri Abhyankar       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr);
175486d161a7SShri Abhyankar     }
175586d161a7SShri Abhyankar     ierr = PetscFree(procsnz);CHKERRQ(ierr);
175686d161a7SShri Abhyankar   } else {
175786d161a7SShri Abhyankar     /* receive numeric values */
1758854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&vals);CHKERRQ(ierr);
175986d161a7SShri Abhyankar 
176086d161a7SShri Abhyankar     /* receive message of values*/
176186d161a7SShri Abhyankar     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr);
176286d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
176386d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
176486d161a7SShri Abhyankar 
176586d161a7SShri Abhyankar     /* insert into matrix */
176686d161a7SShri Abhyankar     jj      = rstart;
176786d161a7SShri Abhyankar     smycols = mycols;
176886d161a7SShri Abhyankar     svals   = vals;
176986d161a7SShri Abhyankar     for (i=0; i<m; i++) {
177086d161a7SShri Abhyankar       ierr     = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
177186d161a7SShri Abhyankar       smycols += ourlens[i];
177286d161a7SShri Abhyankar       svals   += ourlens[i];
177386d161a7SShri Abhyankar       jj++;
177486d161a7SShri Abhyankar     }
177586d161a7SShri Abhyankar   }
177649467e7dSBarry Smith   ierr = PetscFree(ourlens);CHKERRQ(ierr);
177786d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
177886d161a7SShri Abhyankar   ierr = PetscFree(mycols);CHKERRQ(ierr);
177986d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
178086d161a7SShri Abhyankar 
178186d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
178286d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
178386d161a7SShri Abhyankar   PetscFunctionReturn(0);
178486d161a7SShri Abhyankar }
178586d161a7SShri Abhyankar 
1786ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
17876e4ee0c6SHong Zhang {
17886e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
17896e4ee0c6SHong Zhang   Mat            a,b;
1790ace3abfcSBarry Smith   PetscBool      flg;
17916e4ee0c6SHong Zhang   PetscErrorCode ierr;
179290ace30eSBarry Smith 
17936e4ee0c6SHong Zhang   PetscFunctionBegin;
17946e4ee0c6SHong Zhang   a    = matA->A;
17956e4ee0c6SHong Zhang   b    = matB->A;
17966e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
1797b2566f29SBarry Smith   ierr = MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
17986e4ee0c6SHong Zhang   PetscFunctionReturn(0);
17996e4ee0c6SHong Zhang }
180090ace30eSBarry Smith 
1801baa3c1c6SHong Zhang PetscErrorCode MatDestroy_MatTransMatMult_MPIDense_MPIDense(Mat A)
1802baa3c1c6SHong Zhang {
1803baa3c1c6SHong Zhang   PetscErrorCode        ierr;
1804baa3c1c6SHong Zhang   Mat_MPIDense          *a = (Mat_MPIDense*)A->data;
1805baa3c1c6SHong Zhang   Mat_TransMatMultDense *atb = a->atbdense;
1806baa3c1c6SHong Zhang 
1807baa3c1c6SHong Zhang   PetscFunctionBegin;
1808c5ef1628SHong Zhang   ierr = PetscFree3(atb->sendbuf,atb->atbarray,atb->recvcounts);CHKERRQ(ierr);
1809baa3c1c6SHong Zhang   ierr = (atb->destroy)(A);CHKERRQ(ierr);
1810baa3c1c6SHong Zhang   ierr = PetscFree(atb);CHKERRQ(ierr);
1811baa3c1c6SHong Zhang   PetscFunctionReturn(0);
1812baa3c1c6SHong Zhang }
1813baa3c1c6SHong Zhang 
1814cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C)
1815cb20be35SHong Zhang {
1816baa3c1c6SHong Zhang   Mat_MPIDense   *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data;
1817660d5466SHong Zhang   Mat_SeqDense   *aseq=(Mat_SeqDense*)(a->A)->data, *bseq=(Mat_SeqDense*)(b->A)->data;
1818baa3c1c6SHong Zhang   Mat_TransMatMultDense *atb = c->atbdense;
1819cb20be35SHong Zhang   PetscErrorCode ierr;
1820cb20be35SHong Zhang   MPI_Comm       comm;
1821d5017740SHong Zhang   PetscMPIInt    rank,size,*recvcounts=atb->recvcounts;
1822c5ef1628SHong Zhang   PetscScalar    *carray,*atbarray=atb->atbarray,*sendbuf=atb->sendbuf;
1823d5017740SHong Zhang   PetscInt       i,cN=C->cmap->N,cM=C->rmap->N,proc,k,j;
1824660d5466SHong Zhang   PetscScalar    _DOne=1.0,_DZero=0.0;
1825adc7a786SHong Zhang   PetscBLASInt   am,an,bn,aN;
1826e68c0b26SHong Zhang   const PetscInt *ranges;
1827cb20be35SHong Zhang 
1828cb20be35SHong Zhang   PetscFunctionBegin;
1829cb20be35SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1830cb20be35SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1831cb20be35SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1832e68c0b26SHong Zhang 
1833c5ef1628SHong Zhang   /* compute atbarray = aseq^T * bseq */
1834660d5466SHong Zhang   ierr = PetscBLASIntCast(a->A->cmap->n,&an);CHKERRQ(ierr);
1835660d5466SHong Zhang   ierr = PetscBLASIntCast(b->A->cmap->n,&bn);CHKERRQ(ierr);
1836660d5466SHong Zhang   ierr = PetscBLASIntCast(a->A->rmap->n,&am);CHKERRQ(ierr);
1837adc7a786SHong Zhang   ierr = PetscBLASIntCast(A->cmap->N,&aN);CHKERRQ(ierr);
1838adc7a786SHong Zhang   PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&an,&bn,&am,&_DOne,aseq->v,&aseq->lda,bseq->v,&bseq->lda,&_DZero,atbarray,&aN));
1839cb20be35SHong Zhang 
1840cb20be35SHong Zhang   ierr = MatGetOwnershipRanges(C,&ranges);CHKERRQ(ierr);
1841c5ef1628SHong Zhang   for (i=0; i<size; i++) recvcounts[i] = (ranges[i+1] - ranges[i])*cN;
1842cb20be35SHong Zhang 
1843660d5466SHong Zhang   /* arrange atbarray into sendbuf */
1844cb20be35SHong Zhang   k = 0;
1845cb20be35SHong Zhang   for (proc=0; proc<size; proc++) {
1846baa3c1c6SHong Zhang     for (j=0; j<cN; j++) {
1847c5ef1628SHong Zhang       for (i=ranges[proc]; i<ranges[proc+1]; i++) sendbuf[k++] = atbarray[i+j*cM];
1848cb20be35SHong Zhang     }
1849cb20be35SHong Zhang   }
1850c5ef1628SHong Zhang   /* sum all atbarray to local values of C */
1851660d5466SHong Zhang   ierr = MatDenseGetArray(c->A,&carray);CHKERRQ(ierr);
18523462b7efSHong Zhang   ierr = MPI_Reduce_scatter(sendbuf,carray,recvcounts,MPIU_SCALAR,MPIU_SUM,comm);CHKERRQ(ierr);
1853660d5466SHong Zhang   ierr = MatDenseRestoreArray(c->A,&carray);CHKERRQ(ierr);
1854cb20be35SHong Zhang   PetscFunctionReturn(0);
1855cb20be35SHong Zhang }
1856cb20be35SHong Zhang 
1857cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
1858cb20be35SHong Zhang {
1859cb20be35SHong Zhang   PetscErrorCode        ierr;
1860baa3c1c6SHong Zhang   Mat                   Cdense;
1861cb20be35SHong Zhang   MPI_Comm              comm;
1862baa3c1c6SHong Zhang   PetscMPIInt           size;
1863660d5466SHong Zhang   PetscInt              cm=A->cmap->n,cM,cN=B->cmap->N;
1864baa3c1c6SHong Zhang   Mat_MPIDense          *c;
1865baa3c1c6SHong Zhang   Mat_TransMatMultDense *atb;
1866cb20be35SHong Zhang 
1867cb20be35SHong Zhang   PetscFunctionBegin;
1868baa3c1c6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1869cb20be35SHong Zhang   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) {
1870cb20be35SHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
1871cb20be35SHong Zhang   }
1872cb20be35SHong Zhang 
1873baa3c1c6SHong Zhang   /* create matrix product Cdense */
1874baa3c1c6SHong Zhang   ierr = MatCreate(comm,&Cdense);CHKERRQ(ierr);
1875660d5466SHong Zhang   ierr = MatSetSizes(Cdense,cm,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1876baa3c1c6SHong Zhang   ierr = MatSetType(Cdense,MATMPIDENSE);CHKERRQ(ierr);
1877baa3c1c6SHong Zhang   ierr = MatMPIDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr);
1878baa3c1c6SHong Zhang   ierr = MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1879baa3c1c6SHong Zhang   ierr = MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1880baa3c1c6SHong Zhang   *C   = Cdense;
1881baa3c1c6SHong Zhang 
1882baa3c1c6SHong Zhang   /* create data structure for reuse Cdense */
1883baa3c1c6SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1884baa3c1c6SHong Zhang   ierr = PetscNew(&atb);CHKERRQ(ierr);
1885660d5466SHong Zhang   cM = Cdense->rmap->N;
1886c5ef1628SHong Zhang   ierr = PetscMalloc3(cM*cN,&atb->sendbuf,cM*cN,&atb->atbarray,size,&atb->recvcounts);CHKERRQ(ierr);
1887baa3c1c6SHong Zhang 
1888baa3c1c6SHong Zhang   c                    = (Mat_MPIDense*)Cdense->data;
1889baa3c1c6SHong Zhang   c->atbdense          = atb;
1890baa3c1c6SHong Zhang   atb->destroy         = Cdense->ops->destroy;
1891baa3c1c6SHong Zhang   Cdense->ops->destroy = MatDestroy_MatTransMatMult_MPIDense_MPIDense;
1892cb20be35SHong Zhang   PetscFunctionReturn(0);
1893cb20be35SHong Zhang }
1894cb20be35SHong Zhang 
1895cb20be35SHong Zhang PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1896cb20be35SHong Zhang {
1897cb20be35SHong Zhang   PetscErrorCode ierr;
1898cb20be35SHong Zhang 
1899cb20be35SHong Zhang   PetscFunctionBegin;
1900cb20be35SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
1901cb20be35SHong Zhang     ierr = MatTransposeMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr);
1902cb20be35SHong Zhang   }
1903cb20be35SHong Zhang   ierr = MatTransposeMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr);
1904cb20be35SHong Zhang   PetscFunctionReturn(0);
1905cb20be35SHong Zhang }
1906320f2790SHong Zhang 
1907320f2790SHong Zhang PetscErrorCode MatDestroy_MatMatMult_MPIDense_MPIDense(Mat A)
1908320f2790SHong Zhang {
1909320f2790SHong Zhang   PetscErrorCode   ierr;
1910320f2790SHong Zhang   Mat_MPIDense     *a = (Mat_MPIDense*)A->data;
1911320f2790SHong Zhang   Mat_MatMultDense *ab = a->abdense;
1912320f2790SHong Zhang 
1913320f2790SHong Zhang   PetscFunctionBegin;
1914320f2790SHong Zhang   ierr = MatDestroy(&ab->Ce);CHKERRQ(ierr);
1915320f2790SHong Zhang   ierr = MatDestroy(&ab->Ae);CHKERRQ(ierr);
1916320f2790SHong Zhang   ierr = MatDestroy(&ab->Be);CHKERRQ(ierr);
1917320f2790SHong Zhang 
1918320f2790SHong Zhang   ierr = (ab->destroy)(A);CHKERRQ(ierr);
1919320f2790SHong Zhang   ierr = PetscFree(ab);CHKERRQ(ierr);
1920320f2790SHong Zhang   PetscFunctionReturn(0);
1921320f2790SHong Zhang }
1922320f2790SHong Zhang 
19235242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
1924320f2790SHong Zhang PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C)
1925320f2790SHong Zhang {
1926320f2790SHong Zhang   PetscErrorCode   ierr;
1927320f2790SHong Zhang   Mat_MPIDense     *c=(Mat_MPIDense*)C->data;
1928320f2790SHong Zhang   Mat_MatMultDense *ab=c->abdense;
1929320f2790SHong Zhang 
1930320f2790SHong Zhang   PetscFunctionBegin;
1931de0a22f0SHong Zhang   ierr = MatConvert_MPIDense_Elemental(A,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Ae);CHKERRQ(ierr);
1932de0a22f0SHong Zhang   ierr = MatConvert_MPIDense_Elemental(B,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Be);CHKERRQ(ierr);
1933320f2790SHong Zhang   ierr = MatMatMultNumeric(ab->Ae,ab->Be,ab->Ce);CHKERRQ(ierr);
1934de0a22f0SHong Zhang   ierr = MatConvert(ab->Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
1935320f2790SHong Zhang   PetscFunctionReturn(0);
1936320f2790SHong Zhang }
1937320f2790SHong Zhang 
1938320f2790SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
1939320f2790SHong Zhang {
1940320f2790SHong Zhang   PetscErrorCode   ierr;
1941320f2790SHong Zhang   Mat              Ae,Be,Ce;
1942320f2790SHong Zhang   Mat_MPIDense     *c;
1943320f2790SHong Zhang   Mat_MatMultDense *ab;
1944320f2790SHong Zhang 
1945320f2790SHong Zhang   PetscFunctionBegin;
1946320f2790SHong Zhang   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) {
1947378336b6SHong Zhang     SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
1948320f2790SHong Zhang   }
1949320f2790SHong Zhang 
1950320f2790SHong Zhang   /* convert A and B to Elemental matrices Ae and Be */
1951320f2790SHong Zhang   ierr = MatConvert(A,MATELEMENTAL,MAT_INITIAL_MATRIX, &Ae);CHKERRQ(ierr);
1952320f2790SHong Zhang   ierr = MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX, &Be);CHKERRQ(ierr);
1953320f2790SHong Zhang 
1954320f2790SHong Zhang   /* Ce = Ae*Be */
1955320f2790SHong Zhang   ierr = MatMatMultSymbolic(Ae,Be,fill,&Ce);CHKERRQ(ierr);
1956320f2790SHong Zhang   ierr = MatMatMultNumeric(Ae,Be,Ce);CHKERRQ(ierr);
1957320f2790SHong Zhang 
1958320f2790SHong Zhang   /* convert Ce to C */
1959320f2790SHong Zhang   ierr = MatConvert(Ce,MATMPIDENSE,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
1960320f2790SHong Zhang 
1961320f2790SHong Zhang   /* create data structure for reuse Cdense */
1962320f2790SHong Zhang   ierr = PetscNew(&ab);CHKERRQ(ierr);
1963320f2790SHong Zhang   c                  = (Mat_MPIDense*)(*C)->data;
1964320f2790SHong Zhang   c->abdense         = ab;
1965320f2790SHong Zhang 
1966320f2790SHong Zhang   ab->Ae             = Ae;
1967320f2790SHong Zhang   ab->Be             = Be;
1968320f2790SHong Zhang   ab->Ce             = Ce;
1969320f2790SHong Zhang   ab->destroy        = (*C)->ops->destroy;
1970320f2790SHong Zhang   (*C)->ops->destroy        = MatDestroy_MatMatMult_MPIDense_MPIDense;
19719775878aSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIDense_MPIDense;
1972320f2790SHong Zhang   PetscFunctionReturn(0);
1973320f2790SHong Zhang }
1974320f2790SHong Zhang 
1975150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1976320f2790SHong Zhang {
1977320f2790SHong Zhang   PetscErrorCode ierr;
1978320f2790SHong Zhang 
1979320f2790SHong Zhang   PetscFunctionBegin;
1980320f2790SHong Zhang   if (scall == MAT_INITIAL_MATRIX) { /* simbolic product includes numeric product */
198157299f2fSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
1982320f2790SHong Zhang     ierr = MatMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr);
198357299f2fSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
1984320f2790SHong Zhang   } else {
198557299f2fSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
1986320f2790SHong Zhang     ierr = MatMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr);
198757299f2fSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
1988320f2790SHong Zhang   }
1989320f2790SHong Zhang   PetscFunctionReturn(0);
1990320f2790SHong Zhang }
19915242a7b1SHong Zhang #endif
1992