xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 8c778c5567ea6c072b54bc0633dc5ef64cbedb91)
1be1d678aSKris Buschelman 
2ed3cc1f0SBarry Smith /*
3ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
4ed3cc1f0SBarry Smith */
5ed3cc1f0SBarry Smith 
604fea9ffSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>    /*I   "petscmat.h"  I*/
88965ea79SLois Curfman McInnes 
9ba8c8a56SBarry Smith #undef __FUNCT__
10ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
11ab92ecdeSBarry Smith /*@
12ab92ecdeSBarry Smith 
13ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
14ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
15ab92ecdeSBarry Smith 
16ab92ecdeSBarry Smith     Input Parameter:
17ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
18ab92ecdeSBarry Smith 
19ab92ecdeSBarry Smith     Output Parameter:
20ab92ecdeSBarry Smith .      B - the inner matrix
21ab92ecdeSBarry Smith 
228e6c10adSSatish Balay     Level: intermediate
238e6c10adSSatish Balay 
24ab92ecdeSBarry Smith @*/
25ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
26ab92ecdeSBarry Smith {
27ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
28ab92ecdeSBarry Smith   PetscErrorCode ierr;
29ace3abfcSBarry Smith   PetscBool      flg;
30ab92ecdeSBarry Smith 
31ab92ecdeSBarry Smith   PetscFunctionBegin;
32251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
33ab92ecdeSBarry Smith   if (flg) {
34ab92ecdeSBarry Smith     *B = mat->A;
35ab92ecdeSBarry Smith   } else {
36ab92ecdeSBarry Smith     *B = A;
37ab92ecdeSBarry Smith   }
38ab92ecdeSBarry Smith   PetscFunctionReturn(0);
39ab92ecdeSBarry Smith }
40ab92ecdeSBarry Smith 
41ab92ecdeSBarry Smith #undef __FUNCT__
42ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
43ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
44ba8c8a56SBarry Smith {
45ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
46ba8c8a56SBarry Smith   PetscErrorCode ierr;
47d0f46423SBarry Smith   PetscInt       lrow,rstart = A->rmap->rstart,rend = A->rmap->rend;
48ba8c8a56SBarry Smith 
49ba8c8a56SBarry Smith   PetscFunctionBegin;
50e7e72b3dSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows");
51ba8c8a56SBarry Smith   lrow = row - rstart;
52ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
53ba8c8a56SBarry Smith   PetscFunctionReturn(0);
54ba8c8a56SBarry Smith }
55ba8c8a56SBarry Smith 
56ba8c8a56SBarry Smith #undef __FUNCT__
57ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
58ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
59ba8c8a56SBarry Smith {
60ba8c8a56SBarry Smith   PetscErrorCode ierr;
61ba8c8a56SBarry Smith 
62ba8c8a56SBarry Smith   PetscFunctionBegin;
63ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
64ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
65ba8c8a56SBarry Smith   PetscFunctionReturn(0);
66ba8c8a56SBarry Smith }
67ba8c8a56SBarry Smith 
680de54da6SSatish Balay EXTERN_C_BEGIN
694a2ae208SSatish Balay #undef __FUNCT__
704a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
7111bd1e4dSLisandro Dalcin PetscErrorCode  MatGetDiagonalBlock_MPIDense(Mat A,Mat *a)
720de54da6SSatish Balay {
730de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
746849ba73SBarry Smith   PetscErrorCode ierr;
75d0f46423SBarry Smith   PetscInt       m = A->rmap->n,rstart = A->rmap->rstart;
7687828ca2SBarry Smith   PetscScalar    *array;
770de54da6SSatish Balay   MPI_Comm       comm;
7811bd1e4dSLisandro Dalcin   Mat            B;
790de54da6SSatish Balay 
800de54da6SSatish Balay   PetscFunctionBegin;
81e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only square matrices supported.");
820de54da6SSatish Balay 
8311bd1e4dSLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);CHKERRQ(ierr);
8411bd1e4dSLisandro Dalcin   if (!B) {
850de54da6SSatish Balay     ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
8611bd1e4dSLisandro Dalcin     ierr = MatCreate(comm,&B);CHKERRQ(ierr);
8711bd1e4dSLisandro Dalcin     ierr = MatSetSizes(B,m,m,m,m);CHKERRQ(ierr);
8811bd1e4dSLisandro Dalcin     ierr = MatSetType(B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr);
89*8c778c55SBarry Smith     ierr = MatDenseGetArray(mdn->A,&array);CHKERRQ(ierr);
9011bd1e4dSLisandro Dalcin     ierr = MatSeqDenseSetPreallocation(B,array+m*rstart);CHKERRQ(ierr);
91*8c778c55SBarry Smith     ierr = MatDenseRestoreArray(mdn->A,&array);CHKERRQ(ierr);
9211bd1e4dSLisandro Dalcin     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9311bd1e4dSLisandro Dalcin     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9411bd1e4dSLisandro Dalcin     ierr = PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);CHKERRQ(ierr);
9511bd1e4dSLisandro Dalcin     *a = B;
9611bd1e4dSLisandro Dalcin     ierr = MatDestroy(&B);CHKERRQ(ierr);
9711bd1e4dSLisandro Dalcin   } else {
9811bd1e4dSLisandro Dalcin     *a = B;
9911bd1e4dSLisandro Dalcin   }
1000de54da6SSatish Balay   PetscFunctionReturn(0);
1010de54da6SSatish Balay }
1020de54da6SSatish Balay EXTERN_C_END
1030de54da6SSatish Balay 
1044a2ae208SSatish Balay #undef __FUNCT__
1054a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
10613f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1078965ea79SLois Curfman McInnes {
10839b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
109dfbe8321SBarry Smith   PetscErrorCode ierr;
110d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
111ace3abfcSBarry Smith   PetscBool      roworiented = A->roworiented;
1128965ea79SLois Curfman McInnes 
1133a40ed3dSBarry Smith   PetscFunctionBegin;
11471fd2e92SBarry Smith   if (v) PetscValidScalarPointer(v,6);
1158965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1165ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
117e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1188965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1198965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
12039b7565bSBarry Smith       if (roworiented) {
12139b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1223a40ed3dSBarry Smith       } else {
1238965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1245ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
125e32f2f54SBarry Smith           if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
12639b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
12739b7565bSBarry Smith         }
1288965ea79SLois Curfman McInnes       }
1293a40ed3dSBarry Smith     } else {
1303782ba37SSatish Balay       if (!A->donotstash) {
1315080c13bSMatthew G Knepley         mat->assembled = PETSC_FALSE;
13239b7565bSBarry Smith         if (roworiented) {
133b400d20cSBarry Smith           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n,PETSC_FALSE);CHKERRQ(ierr);
134d36fbae8SSatish Balay         } else {
135b400d20cSBarry Smith           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m,PETSC_FALSE);CHKERRQ(ierr);
13639b7565bSBarry Smith         }
137b49de8d1SLois Curfman McInnes       }
138b49de8d1SLois Curfman McInnes     }
1393782ba37SSatish Balay   }
1403a40ed3dSBarry Smith   PetscFunctionReturn(0);
141b49de8d1SLois Curfman McInnes }
142b49de8d1SLois Curfman McInnes 
1434a2ae208SSatish Balay #undef __FUNCT__
1444a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
14513f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
146b49de8d1SLois Curfman McInnes {
147b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
148dfbe8321SBarry Smith   PetscErrorCode ierr;
149d0f46423SBarry Smith   PetscInt       i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row;
150b49de8d1SLois Curfman McInnes 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
152b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
153e32f2f54SBarry Smith     if (idxm[i] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
154e32f2f54SBarry Smith     if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
155b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
156b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
157b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
158e32f2f54SBarry Smith         if (idxn[j] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
159e7e72b3dSBarry Smith         if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
160b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
161b49de8d1SLois Curfman McInnes       }
162e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
1638965ea79SLois Curfman McInnes   }
1643a40ed3dSBarry Smith   PetscFunctionReturn(0);
1658965ea79SLois Curfman McInnes }
1668965ea79SLois Curfman McInnes 
1674a2ae208SSatish Balay #undef __FUNCT__
168*8c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_MPIDense"
169*8c778c55SBarry Smith PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[])
170ff14e315SSatish Balay {
171ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
172dfbe8321SBarry Smith   PetscErrorCode ierr;
173ff14e315SSatish Balay 
1743a40ed3dSBarry Smith   PetscFunctionBegin;
175*8c778c55SBarry Smith   ierr = MatDenseGetArray(a->A,array);CHKERRQ(ierr);
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177ff14e315SSatish Balay }
178ff14e315SSatish Balay 
1794a2ae208SSatish Balay #undef __FUNCT__
1804a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
1814aa3045dSJed Brown static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
182ca3fa75bSLois Curfman McInnes {
183ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data,*newmatd;
184ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
1856849ba73SBarry Smith   PetscErrorCode ierr;
1864aa3045dSJed Brown   PetscInt       i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols;
1875d0c19d7SBarry Smith   const PetscInt *irow,*icol;
18887828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
189ca3fa75bSLois Curfman McInnes   Mat            newmat;
1904aa3045dSJed Brown   IS             iscol_local;
191ca3fa75bSLois Curfman McInnes 
192ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
1934aa3045dSJed Brown   ierr = ISAllGather(iscol,&iscol_local);CHKERRQ(ierr);
194ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
1954aa3045dSJed Brown   ierr = ISGetIndices(iscol_local,&icol);CHKERRQ(ierr);
196b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
197b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
1984aa3045dSJed Brown   ierr = ISGetSize(iscol,&Ncols);CHKERRQ(ierr); /* global number of columns, size of iscol_local */
199ca3fa75bSLois Curfman McInnes 
200ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
2017eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
2027eba5e9cSLois Curfman McInnes      original matrix! */
203ca3fa75bSLois Curfman McInnes 
204ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
2057eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
206ca3fa75bSLois Curfman McInnes 
207ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
208ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
209e32f2f54SBarry Smith     /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2107eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
211ca3fa75bSLois Curfman McInnes     newmat = *B;
212ca3fa75bSLois Curfman McInnes   } else {
213ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
2147adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
2154aa3045dSJed Brown     ierr = MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);CHKERRQ(ierr);
2167adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
217878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
218ca3fa75bSLois Curfman McInnes   }
219ca3fa75bSLois Curfman McInnes 
220ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
221ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
222ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
223ca3fa75bSLois Curfman McInnes 
2244aa3045dSJed Brown   for (i=0; i<Ncols; i++) {
22525a33276SHong Zhang     av = v + ((Mat_SeqDense *)mat->A->data)->lda*icol[i];
226ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2277eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
228ca3fa75bSLois Curfman McInnes     }
229ca3fa75bSLois Curfman McInnes   }
230ca3fa75bSLois Curfman McInnes 
231ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
232ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
233ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
234ca3fa75bSLois Curfman McInnes 
235ca3fa75bSLois Curfman McInnes   /* Free work space */
236ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
2375bdf786aSShri Abhyankar   ierr = ISRestoreIndices(iscol_local,&icol);CHKERRQ(ierr);
23832bb1f2dSLisandro Dalcin   ierr = ISDestroy(&iscol_local);CHKERRQ(ierr);
239ca3fa75bSLois Curfman McInnes   *B = newmat;
240ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
241ca3fa75bSLois Curfman McInnes }
242ca3fa75bSLois Curfman McInnes 
2434a2ae208SSatish Balay #undef __FUNCT__
244*8c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_MPIDense"
245*8c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[])
246ff14e315SSatish Balay {
24773a71a0fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
24873a71a0fSBarry Smith   PetscErrorCode ierr;
24973a71a0fSBarry Smith 
2503a40ed3dSBarry Smith   PetscFunctionBegin;
251*8c778c55SBarry Smith   ierr = MatDenseRestoreArray(a->A,array);CHKERRQ(ierr);
2523a40ed3dSBarry Smith   PetscFunctionReturn(0);
253ff14e315SSatish Balay }
254ff14e315SSatish Balay 
2554a2ae208SSatish Balay #undef __FUNCT__
2564a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
257dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2588965ea79SLois Curfman McInnes {
25939ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
2607adad957SLisandro Dalcin   MPI_Comm       comm = ((PetscObject)mat)->comm;
261dfbe8321SBarry Smith   PetscErrorCode ierr;
26213f74950SBarry Smith   PetscInt       nstash,reallocs;
2638965ea79SLois Curfman McInnes   InsertMode     addv;
2648965ea79SLois Curfman McInnes 
2653a40ed3dSBarry Smith   PetscFunctionBegin;
2668965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
267ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
268e7e72b3dSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
269e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2708965ea79SLois Curfman McInnes 
271d0f46423SBarry Smith   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr);
2728798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
273ae15b995SBarry Smith   ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
2743a40ed3dSBarry Smith   PetscFunctionReturn(0);
2758965ea79SLois Curfman McInnes }
2768965ea79SLois Curfman McInnes 
2774a2ae208SSatish Balay #undef __FUNCT__
2784a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
279dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2808965ea79SLois Curfman McInnes {
28139ddd567SLois Curfman McInnes   Mat_MPIDense    *mdn=(Mat_MPIDense*)mat->data;
2826849ba73SBarry Smith   PetscErrorCode  ierr;
28313f74950SBarry Smith   PetscInt        i,*row,*col,flg,j,rstart,ncols;
28413f74950SBarry Smith   PetscMPIInt     n;
28587828ca2SBarry Smith   PetscScalar     *val;
286e0fa3b82SLois Curfman McInnes   InsertMode      addv=mat->insertmode;
2878965ea79SLois Curfman McInnes 
2883a40ed3dSBarry Smith   PetscFunctionBegin;
2898965ea79SLois Curfman McInnes   /*  wait on receives */
2907ef1d9bdSSatish Balay   while (1) {
2918798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2927ef1d9bdSSatish Balay     if (!flg) break;
2938965ea79SLois Curfman McInnes 
2947ef1d9bdSSatish Balay     for (i=0; i<n;) {
2957ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2967ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2977ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2987ef1d9bdSSatish Balay       else       ncols = n-i;
2997ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
3007ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
3017ef1d9bdSSatish Balay       i = j;
3028965ea79SLois Curfman McInnes     }
3037ef1d9bdSSatish Balay   }
3048798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
3058965ea79SLois Curfman McInnes 
30639ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
30739ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3088965ea79SLois Curfman McInnes 
3096d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
31039ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3118965ea79SLois Curfman McInnes   }
3123a40ed3dSBarry Smith   PetscFunctionReturn(0);
3138965ea79SLois Curfman McInnes }
3148965ea79SLois Curfman McInnes 
3154a2ae208SSatish Balay #undef __FUNCT__
3164a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
317dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3188965ea79SLois Curfman McInnes {
319dfbe8321SBarry Smith   PetscErrorCode ierr;
32039ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3213a40ed3dSBarry Smith 
3223a40ed3dSBarry Smith   PetscFunctionBegin;
3233a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3243a40ed3dSBarry Smith   PetscFunctionReturn(0);
3258965ea79SLois Curfman McInnes }
3268965ea79SLois Curfman McInnes 
3278965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3288965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3298965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3303501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3318965ea79SLois Curfman McInnes    routine.
3328965ea79SLois Curfman McInnes */
3334a2ae208SSatish Balay #undef __FUNCT__
3344a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
3352b40b63fSBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
3368965ea79SLois Curfman McInnes {
33739ddd567SLois Curfman McInnes   Mat_MPIDense      *l = (Mat_MPIDense*)A->data;
3386849ba73SBarry Smith   PetscErrorCode    ierr;
339d0f46423SBarry Smith   PetscInt          i,*owners = A->rmap->range;
34013f74950SBarry Smith   PetscInt          *nprocs,j,idx,nsends;
34113f74950SBarry Smith   PetscInt          nmax,*svalues,*starts,*owner,nrecvs;
3427adad957SLisandro Dalcin   PetscInt          *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source;
34313f74950SBarry Smith   PetscInt          *lens,*lrows,*values;
34413f74950SBarry Smith   PetscMPIInt       n,imdex,rank = l->rank,size = l->size;
3457adad957SLisandro Dalcin   MPI_Comm          comm = ((PetscObject)A)->comm;
3468965ea79SLois Curfman McInnes   MPI_Request       *send_waits,*recv_waits;
3478965ea79SLois Curfman McInnes   MPI_Status        recv_status,*send_status;
348ace3abfcSBarry Smith   PetscBool         found;
34997b48c8fSBarry Smith   const PetscScalar *xx;
35097b48c8fSBarry Smith   PetscScalar       *bb;
3518965ea79SLois Curfman McInnes 
3523a40ed3dSBarry Smith   PetscFunctionBegin;
353b9679d65SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Only handles square matrices");
354b9679d65SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership");
3558965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
35613f74950SBarry Smith   ierr  = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
35713f74950SBarry Smith   ierr  = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
35813f74950SBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/
3598965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3608965ea79SLois Curfman McInnes     idx = rows[i];
36135d8aa7fSBarry Smith     found = PETSC_FALSE;
3628965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3638965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
364c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3658965ea79SLois Curfman McInnes       }
3668965ea79SLois Curfman McInnes     }
367e32f2f54SBarry Smith     if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3688965ea79SLois Curfman McInnes   }
369c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3708965ea79SLois Curfman McInnes 
3718965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
372c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3738965ea79SLois Curfman McInnes 
3748965ea79SLois Curfman McInnes   /* post receives:   */
37513f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
376b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
37813f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3798965ea79SLois Curfman McInnes   }
3808965ea79SLois Curfman McInnes 
3818965ea79SLois Curfman McInnes   /* do sends:
3828965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3838965ea79SLois Curfman McInnes          the ith processor
3848965ea79SLois Curfman McInnes   */
38513f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
386b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
38713f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
3888965ea79SLois Curfman McInnes   starts[0]  = 0;
389c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3908965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3918965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3928965ea79SLois Curfman McInnes   }
3938965ea79SLois Curfman McInnes 
3948965ea79SLois Curfman McInnes   starts[0] = 0;
395c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3968965ea79SLois Curfman McInnes   count = 0;
3978965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
398c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
39913f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
4008965ea79SLois Curfman McInnes     }
4018965ea79SLois Curfman McInnes   }
402606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
4038965ea79SLois Curfman McInnes 
4048965ea79SLois Curfman McInnes   base = owners[rank];
4058965ea79SLois Curfman McInnes 
4068965ea79SLois Curfman McInnes   /*  wait on receives */
40774ed9c26SBarry Smith   ierr   = PetscMalloc2(nrecvs,PetscInt,&lens,nrecvs,PetscInt,&source);CHKERRQ(ierr);
40874ed9c26SBarry Smith   count  = nrecvs;
40974ed9c26SBarry Smith   slen   = 0;
4108965ea79SLois Curfman McInnes   while (count) {
411ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4128965ea79SLois Curfman McInnes     /* unpack receives into our local space */
41313f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4148965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
4158965ea79SLois Curfman McInnes     lens[imdex]    = n;
4168965ea79SLois Curfman McInnes     slen += n;
4178965ea79SLois Curfman McInnes     count--;
4188965ea79SLois Curfman McInnes   }
419606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4208965ea79SLois Curfman McInnes 
4218965ea79SLois Curfman McInnes   /* move the data into the send scatter */
42213f74950SBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
4238965ea79SLois Curfman McInnes   count = 0;
4248965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4258965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4268965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4278965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4288965ea79SLois Curfman McInnes     }
4298965ea79SLois Curfman McInnes   }
430606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
43174ed9c26SBarry Smith   ierr = PetscFree2(lens,source);CHKERRQ(ierr);
432606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
433606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
4348965ea79SLois Curfman McInnes 
43597b48c8fSBarry Smith   /* fix right hand side if needed */
43697b48c8fSBarry Smith   if (x && b) {
43797b48c8fSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
43897b48c8fSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
43997b48c8fSBarry Smith     for (i=0; i<slen; i++) {
44097b48c8fSBarry Smith       bb[lrows[i]] = diag*xx[lrows[i]];
44197b48c8fSBarry Smith     }
44297b48c8fSBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
44397b48c8fSBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
44497b48c8fSBarry Smith   }
44597b48c8fSBarry Smith 
4468965ea79SLois Curfman McInnes   /* actually zap the local rows */
447b9679d65SBarry Smith   ierr = MatZeroRows(l->A,slen,lrows,0.0,0,0);CHKERRQ(ierr);
448e2eb51b1SBarry Smith   if (diag != 0.0) {
449b9679d65SBarry Smith     Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data;
450b9679d65SBarry Smith     PetscInt      m = ll->lda, i;
451b9679d65SBarry Smith 
452b9679d65SBarry Smith     for (i=0; i<slen; i++) {
453b9679d65SBarry Smith       ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag;
454b9679d65SBarry Smith     }
455b9679d65SBarry Smith   }
456606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4578965ea79SLois Curfman McInnes 
4588965ea79SLois Curfman McInnes   /* wait on sends */
4598965ea79SLois Curfman McInnes   if (nsends) {
460b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
461ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
462606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4638965ea79SLois Curfman McInnes   }
464606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
465606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4668965ea79SLois Curfman McInnes 
4673a40ed3dSBarry Smith   PetscFunctionReturn(0);
4688965ea79SLois Curfman McInnes }
4698965ea79SLois Curfman McInnes 
4704a2ae208SSatish Balay #undef __FUNCT__
4714a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
472dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4738965ea79SLois Curfman McInnes {
47439ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
475dfbe8321SBarry Smith   PetscErrorCode ierr;
476c456f294SBarry Smith 
4773a40ed3dSBarry Smith   PetscFunctionBegin;
478ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
479ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
48044cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4813a40ed3dSBarry Smith   PetscFunctionReturn(0);
4828965ea79SLois Curfman McInnes }
4838965ea79SLois Curfman McInnes 
4844a2ae208SSatish Balay #undef __FUNCT__
4854a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
486dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4878965ea79SLois Curfman McInnes {
48839ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
489dfbe8321SBarry Smith   PetscErrorCode ierr;
490c456f294SBarry Smith 
4913a40ed3dSBarry Smith   PetscFunctionBegin;
492ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
493ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
49444cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4953a40ed3dSBarry Smith   PetscFunctionReturn(0);
4968965ea79SLois Curfman McInnes }
4978965ea79SLois Curfman McInnes 
4984a2ae208SSatish Balay #undef __FUNCT__
4994a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
500dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
501096963f5SLois Curfman McInnes {
502096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
503dfbe8321SBarry Smith   PetscErrorCode ierr;
50487828ca2SBarry Smith   PetscScalar    zero = 0.0;
505096963f5SLois Curfman McInnes 
5063a40ed3dSBarry Smith   PetscFunctionBegin;
5072dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
5087c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
509ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
510ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5113a40ed3dSBarry Smith   PetscFunctionReturn(0);
512096963f5SLois Curfman McInnes }
513096963f5SLois Curfman McInnes 
5144a2ae208SSatish Balay #undef __FUNCT__
5154a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
516dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
517096963f5SLois Curfman McInnes {
518096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
519dfbe8321SBarry Smith   PetscErrorCode ierr;
520096963f5SLois Curfman McInnes 
5213a40ed3dSBarry Smith   PetscFunctionBegin;
5223501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
5237c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
524ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
525ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5263a40ed3dSBarry Smith   PetscFunctionReturn(0);
527096963f5SLois Curfman McInnes }
528096963f5SLois Curfman McInnes 
5294a2ae208SSatish Balay #undef __FUNCT__
5304a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
531dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5328965ea79SLois Curfman McInnes {
53339ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
534096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
535dfbe8321SBarry Smith   PetscErrorCode ierr;
536d0f46423SBarry Smith   PetscInt       len,i,n,m = A->rmap->n,radd;
53787828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
538ed3cc1f0SBarry Smith 
5393a40ed3dSBarry Smith   PetscFunctionBegin;
5402dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5411ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
542096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
543e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
544d0f46423SBarry Smith   len  = PetscMin(a->A->rmap->n,a->A->cmap->n);
545d0f46423SBarry Smith   radd = A->rmap->rstart*m;
54644cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
547096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
548096963f5SLois Curfman McInnes   }
5491ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5503a40ed3dSBarry Smith   PetscFunctionReturn(0);
5518965ea79SLois Curfman McInnes }
5528965ea79SLois Curfman McInnes 
5534a2ae208SSatish Balay #undef __FUNCT__
5544a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
555dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5568965ea79SLois Curfman McInnes {
5573501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
558dfbe8321SBarry Smith   PetscErrorCode ierr;
559ed3cc1f0SBarry Smith 
5603a40ed3dSBarry Smith   PetscFunctionBegin;
56194d884c6SBarry Smith 
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);
572901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
573901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5744ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5754ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5764ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5773a40ed3dSBarry Smith   PetscFunctionReturn(0);
5788965ea79SLois Curfman McInnes }
57939ddd567SLois Curfman McInnes 
5804a2ae208SSatish Balay #undef __FUNCT__
5814a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5826849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5838965ea79SLois Curfman McInnes {
58439ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
585dfbe8321SBarry Smith   PetscErrorCode    ierr;
586aa05aa95SBarry Smith   PetscViewerFormat format;
587aa05aa95SBarry Smith   int               fd;
588d0f46423SBarry Smith   PetscInt          header[4],mmax,N = mat->cmap->N,i,j,m,k;
589aa05aa95SBarry Smith   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
590578230a0SSatish Balay   PetscScalar       *work,*v,*vv;
591aa05aa95SBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;
5927056b6fcSBarry Smith 
5933a40ed3dSBarry Smith   PetscFunctionBegin;
59439ddd567SLois Curfman McInnes   if (mdn->size == 1) {
59539ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
596aa05aa95SBarry Smith   } else {
597aa05aa95SBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
5987adad957SLisandro Dalcin     ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&rank);CHKERRQ(ierr);
5997adad957SLisandro Dalcin     ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr);
600aa05aa95SBarry Smith 
601aa05aa95SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
602f4403165SShri Abhyankar     if (format == PETSC_VIEWER_NATIVE) {
603aa05aa95SBarry Smith 
604aa05aa95SBarry Smith       if (!rank) {
605aa05aa95SBarry Smith         /* store the matrix as a dense matrix */
6060700a824SBarry Smith         header[0] = MAT_FILE_CLASSID;
607d0f46423SBarry Smith         header[1] = mat->rmap->N;
608aa05aa95SBarry Smith         header[2] = N;
609aa05aa95SBarry Smith         header[3] = MATRIX_BINARY_FORMAT_DENSE;
610aa05aa95SBarry Smith         ierr = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
611aa05aa95SBarry Smith 
612aa05aa95SBarry Smith         /* get largest work array needed for transposing array */
613d0f46423SBarry Smith         mmax = mat->rmap->n;
614aa05aa95SBarry Smith         for (i=1; i<size; i++) {
615d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
6168965ea79SLois Curfman McInnes         }
617aa05aa95SBarry Smith         ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&work);CHKERRQ(ierr);
618aa05aa95SBarry Smith 
619aa05aa95SBarry Smith         /* write out local array, by rows */
620d0f46423SBarry Smith         m    = mat->rmap->n;
621aa05aa95SBarry Smith         v    = a->v;
622aa05aa95SBarry Smith         for (j=0; j<N; j++) {
623aa05aa95SBarry Smith           for (i=0; i<m; i++) {
624578230a0SSatish Balay             work[j + i*N] = *v++;
625aa05aa95SBarry Smith           }
626aa05aa95SBarry Smith         }
627aa05aa95SBarry Smith         ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
628aa05aa95SBarry Smith         /* get largest work array to receive messages from other processes, excludes process zero */
629aa05aa95SBarry Smith         mmax = 0;
630aa05aa95SBarry Smith         for (i=1; i<size; i++) {
631d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
632aa05aa95SBarry Smith         }
633578230a0SSatish Balay         ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&vv);CHKERRQ(ierr);
634aa05aa95SBarry Smith         for(k = 1; k < size; k++) {
635f8009846SMatthew Knepley           v    = vv;
636d0f46423SBarry Smith           m    = mat->rmap->range[k+1] - mat->rmap->range[k];
637a25532f0SBarry Smith           ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,((PetscObject)mat)->comm);CHKERRQ(ierr);
638aa05aa95SBarry Smith 
639aa05aa95SBarry Smith           for(j = 0; j < N; j++) {
640aa05aa95SBarry Smith             for(i = 0; i < m; i++) {
641578230a0SSatish Balay               work[j + i*N] = *v++;
642aa05aa95SBarry Smith             }
643aa05aa95SBarry Smith           }
644aa05aa95SBarry Smith           ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
645aa05aa95SBarry Smith         }
646aa05aa95SBarry Smith         ierr = PetscFree(work);CHKERRQ(ierr);
647578230a0SSatish Balay         ierr = PetscFree(vv);CHKERRQ(ierr);
648aa05aa95SBarry Smith       } else {
649a25532f0SBarry Smith         ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,((PetscObject)mat)->comm);CHKERRQ(ierr);
650aa05aa95SBarry Smith       }
651eb3f98d2SBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE)");
652aa05aa95SBarry Smith   }
6533a40ed3dSBarry Smith   PetscFunctionReturn(0);
6548965ea79SLois Curfman McInnes }
6558965ea79SLois Curfman McInnes 
6564a2ae208SSatish Balay #undef __FUNCT__
6574a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
6586849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6598965ea79SLois Curfman McInnes {
66039ddd567SLois Curfman McInnes   Mat_MPIDense          *mdn = (Mat_MPIDense*)mat->data;
661dfbe8321SBarry Smith   PetscErrorCode        ierr;
66232dcc486SBarry Smith   PetscMPIInt           size = mdn->size,rank = mdn->rank;
663a313700dSBarry Smith   const PetscViewerType vtype;
664ace3abfcSBarry Smith   PetscBool             iascii,isdraw;
665b0a32e0cSBarry Smith   PetscViewer           sviewer;
666f3ef73ceSBarry Smith   PetscViewerFormat     format;
6678965ea79SLois Curfman McInnes 
6683a40ed3dSBarry Smith   PetscFunctionBegin;
669251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
670251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
67132077d6dSBarry Smith   if (iascii) {
672b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
673b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
674456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6754e220ebcSLois Curfman McInnes       MatInfo info;
676888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
6777b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
6787b23a99aSBarry 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);
679b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6807b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
6815d00a290SHong Zhang       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
6823a40ed3dSBarry Smith       PetscFunctionReturn(0);
683fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6843a40ed3dSBarry Smith       PetscFunctionReturn(0);
6858965ea79SLois Curfman McInnes     }
686f1af5d2fSBarry Smith   } else if (isdraw) {
687b0a32e0cSBarry Smith     PetscDraw  draw;
688ace3abfcSBarry Smith     PetscBool  isnull;
689f1af5d2fSBarry Smith 
690b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
691b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
692f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
693f1af5d2fSBarry Smith   }
69477ed5343SBarry Smith 
6958965ea79SLois Curfman McInnes   if (size == 1) {
69639ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6973a40ed3dSBarry Smith   } else {
6988965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
6998965ea79SLois Curfman McInnes     Mat         A;
700d0f46423SBarry Smith     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
701ba8c8a56SBarry Smith     PetscInt    *cols;
702ba8c8a56SBarry Smith     PetscScalar *vals;
7038965ea79SLois Curfman McInnes 
7047adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)mat)->comm,&A);CHKERRQ(ierr);
7058965ea79SLois Curfman McInnes     if (!rank) {
706f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
7073a40ed3dSBarry Smith     } else {
708f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
7098965ea79SLois Curfman McInnes     }
7107adad957SLisandro Dalcin     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
711878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
712878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
71352e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
7148965ea79SLois Curfman McInnes 
71539ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
71639ddd567SLois Curfman McInnes        but it's quick for now */
71751022da4SBarry Smith     A->insertmode = INSERT_VALUES;
718d0f46423SBarry Smith     row = mat->rmap->rstart; m = mdn->A->rmap->n;
7198965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
720ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
721ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
722ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
72339ddd567SLois Curfman McInnes       row++;
7248965ea79SLois Curfman McInnes     }
7258965ea79SLois Curfman McInnes 
7266d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7276d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
728b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
729b9b97703SBarry Smith     if (!rank) {
7307566de4bSShri Abhyankar       ierr = PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr);
7317566de4bSShri Abhyankar       /* Set the type name to MATMPIDense so that the correct type can be printed out by PetscObjectPrintClassNamePrefixType() in MatView_SeqDense_ASCII()*/
7327566de4bSShri Abhyankar       PetscStrcpy(((PetscObject)((Mat_MPIDense*)(A->data))->A)->type_name,MATMPIDENSE);
7336831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
7348965ea79SLois Curfman McInnes     }
735b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
736b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7376bf464f9SBarry Smith     ierr = MatDestroy(&A);CHKERRQ(ierr);
7388965ea79SLois Curfman McInnes   }
7393a40ed3dSBarry Smith   PetscFunctionReturn(0);
7408965ea79SLois Curfman McInnes }
7418965ea79SLois Curfman McInnes 
7424a2ae208SSatish Balay #undef __FUNCT__
7434a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
744dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
7458965ea79SLois Curfman McInnes {
746dfbe8321SBarry Smith   PetscErrorCode ierr;
747ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw,issocket;
7488965ea79SLois Curfman McInnes 
749433994e6SBarry Smith   PetscFunctionBegin;
7500f5bd95cSBarry Smith 
751251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
752251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
753251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr);
754251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
7550f5bd95cSBarry Smith 
75632077d6dSBarry Smith   if (iascii || issocket || isdraw) {
757f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7580f5bd95cSBarry Smith   } else if (isbinary) {
7593a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
760eb3f98d2SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
7613a40ed3dSBarry Smith   PetscFunctionReturn(0);
7628965ea79SLois Curfman McInnes }
7638965ea79SLois Curfman McInnes 
7644a2ae208SSatish Balay #undef __FUNCT__
7654a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
766dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7678965ea79SLois Curfman McInnes {
7683501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7693501a2bdSLois Curfman McInnes   Mat            mdn = mat->A;
770dfbe8321SBarry Smith   PetscErrorCode ierr;
771329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7728965ea79SLois Curfman McInnes 
7733a40ed3dSBarry Smith   PetscFunctionBegin;
7744e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
7754e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7764e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7774e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7788965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7794e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7804e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7814e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7824e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7834e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7848965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
785d9822059SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,((PetscObject)A)->comm);CHKERRQ(ierr);
7864e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7874e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7884e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7894e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7904e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7918965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
792d9822059SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
7934e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7944e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7954e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7964e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7974e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7988965ea79SLois Curfman McInnes   }
7994e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
8004e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
8014e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
8023a40ed3dSBarry Smith   PetscFunctionReturn(0);
8038965ea79SLois Curfman McInnes }
8048965ea79SLois Curfman McInnes 
8054a2ae208SSatish Balay #undef __FUNCT__
8064a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
807ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool  flg)
8088965ea79SLois Curfman McInnes {
80939ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
810dfbe8321SBarry Smith   PetscErrorCode ierr;
8118965ea79SLois Curfman McInnes 
8123a40ed3dSBarry Smith   PetscFunctionBegin;
81312c028f9SKris Buschelman   switch (op) {
814512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
81512c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
81612c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
8174e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
81812c028f9SKris Buschelman     break;
81912c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
8204e0d8c25SBarry Smith     a->roworiented = flg;
8214e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82212c028f9SKris Buschelman     break;
8234e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
82413fa8e87SLisandro Dalcin   case MAT_KEEP_NONZERO_PATTERN:
82512c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
826290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
82712c028f9SKris Buschelman     break;
82812c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
8294e0d8c25SBarry Smith     a->donotstash = flg;
83012c028f9SKris Buschelman     break;
83177e54ba9SKris Buschelman   case MAT_SYMMETRIC:
83277e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8339a4540c5SBarry Smith   case MAT_HERMITIAN:
8349a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
835600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
836290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
83777e54ba9SKris Buschelman     break;
83812c028f9SKris Buschelman   default:
839e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
8403a40ed3dSBarry Smith   }
8413a40ed3dSBarry Smith   PetscFunctionReturn(0);
8428965ea79SLois Curfman McInnes }
8438965ea79SLois Curfman McInnes 
8448965ea79SLois Curfman McInnes 
8454a2ae208SSatish Balay #undef __FUNCT__
8464a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
847dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8485b2fa520SLois Curfman McInnes {
8495b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8505b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
85187828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
852dfbe8321SBarry Smith   PetscErrorCode ierr;
853d0f46423SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n;
8545b2fa520SLois Curfman McInnes 
8555b2fa520SLois Curfman McInnes   PetscFunctionBegin;
85672d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8575b2fa520SLois Curfman McInnes   if (ll) {
85872d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
859e32f2f54SBarry Smith     if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8601ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8615b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8625b2fa520SLois Curfman McInnes       x = l[i];
8635b2fa520SLois Curfman McInnes       v = mat->v + i;
8645b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8655b2fa520SLois Curfman McInnes     }
8661ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
867efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8685b2fa520SLois Curfman McInnes   }
8695b2fa520SLois Curfman McInnes   if (rr) {
870175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
871e32f2f54SBarry Smith     if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
872ca9f406cSSatish Balay     ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873ca9f406cSSatish Balay     ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8741ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8755b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8765b2fa520SLois Curfman McInnes       x = r[i];
8775b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8785b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
8795b2fa520SLois Curfman McInnes     }
8801ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
881efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8825b2fa520SLois Curfman McInnes   }
8835b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8845b2fa520SLois Curfman McInnes }
8855b2fa520SLois Curfman McInnes 
8864a2ae208SSatish Balay #undef __FUNCT__
8874a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
888dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
889096963f5SLois Curfman McInnes {
8903501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8913501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
892dfbe8321SBarry Smith   PetscErrorCode ierr;
89313f74950SBarry Smith   PetscInt       i,j;
894329f5518SBarry Smith   PetscReal      sum = 0.0;
89587828ca2SBarry Smith   PetscScalar    *v = mat->v;
8963501a2bdSLois Curfman McInnes 
8973a40ed3dSBarry Smith   PetscFunctionBegin;
8983501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
899064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
9003501a2bdSLois Curfman McInnes   } else {
9013501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
902d0f46423SBarry Smith       for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) {
903aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
904329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
9053501a2bdSLois Curfman McInnes #else
9063501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
9073501a2bdSLois Curfman McInnes #endif
9083501a2bdSLois Curfman McInnes       }
909d9822059SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,((PetscObject)A)->comm);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;
91474ed9c26SBarry Smith       ierr = PetscMalloc2(A->cmap->N,PetscReal,&tmp,A->cmap->N,PetscReal,&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       }
924d9822059SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
925d0f46423SBarry Smith       for (j=0; j<A->cmap->N; j++) {
926064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
9273501a2bdSLois Curfman McInnes       }
92874ed9c26SBarry Smith       ierr = PetscFree2(tmp,tmp);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);
933d9822059SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,((PetscObject)A)->comm);CHKERRQ(ierr);
934e7e72b3dSBarry Smith     } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support for two norm");
9353501a2bdSLois Curfman McInnes   }
9363a40ed3dSBarry Smith   PetscFunctionReturn(0);
9373501a2bdSLois Curfman McInnes }
9383501a2bdSLois Curfman McInnes 
9394a2ae208SSatish Balay #undef __FUNCT__
9404a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
941fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
9423501a2bdSLois Curfman McInnes {
9433501a2bdSLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
9443501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
9453501a2bdSLois Curfman McInnes   Mat            B;
946d0f46423SBarry Smith   PetscInt       M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart;
9476849ba73SBarry Smith   PetscErrorCode ierr;
94813f74950SBarry Smith   PetscInt       j,i;
94987828ca2SBarry Smith   PetscScalar    *v;
9503501a2bdSLois Curfman McInnes 
9513a40ed3dSBarry Smith   PetscFunctionBegin;
952e7e72b3dSBarry Smith   if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Supports square matrix only in-place");
953fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
9547adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
955d0f46423SBarry Smith     ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr);
9567adad957SLisandro Dalcin     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
957878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
958fc4dec0aSBarry Smith   } else {
959fc4dec0aSBarry Smith     B = *matout;
960fc4dec0aSBarry Smith   }
9613501a2bdSLois Curfman McInnes 
962d0f46423SBarry Smith   m = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v;
9631acff37aSSatish Balay   ierr = PetscMalloc(m*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
9643501a2bdSLois Curfman McInnes   for (i=0; i<m; i++) rwork[i] = rstart + i;
9651acff37aSSatish Balay   for (j=0; j<n; j++) {
9663501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9673501a2bdSLois Curfman McInnes     v   += m;
9683501a2bdSLois Curfman McInnes   }
969606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9706d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9716d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
972815cbec1SBarry Smith   if (reuse == MAT_INITIAL_MATRIX || *matout != A) {
9733501a2bdSLois Curfman McInnes     *matout = B;
9743501a2bdSLois Curfman McInnes   } else {
975eb6b5d47SBarry Smith     ierr = MatHeaderMerge(A,B);CHKERRQ(ierr);
9763501a2bdSLois Curfman McInnes   }
9773a40ed3dSBarry Smith   PetscFunctionReturn(0);
978096963f5SLois Curfman McInnes }
979096963f5SLois Curfman McInnes 
98044cd7ae7SLois Curfman McInnes 
9816849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
982d84338a6SBarry Smith extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar);
9838965ea79SLois Curfman McInnes 
9844a2ae208SSatish Balay #undef __FUNCT__
9854994cf47SJed Brown #define __FUNCT__ "MatSetUp_MPIDense"
9864994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A)
987273d9f13SBarry Smith {
988dfbe8321SBarry Smith   PetscErrorCode ierr;
989273d9f13SBarry Smith 
990273d9f13SBarry Smith   PetscFunctionBegin;
991273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
992273d9f13SBarry Smith   PetscFunctionReturn(0);
993273d9f13SBarry Smith }
994273d9f13SBarry Smith 
99530716080SHong Zhang #undef __FUNCT__
996488007eeSBarry Smith #define __FUNCT__ "MatAXPY_MPIDense"
997488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
998488007eeSBarry Smith {
999488007eeSBarry Smith   PetscErrorCode ierr;
1000488007eeSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;
1001488007eeSBarry Smith 
1002488007eeSBarry Smith   PetscFunctionBegin;
1003488007eeSBarry Smith   ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr);
1004488007eeSBarry Smith   PetscFunctionReturn(0);
1005488007eeSBarry Smith }
1006488007eeSBarry Smith 
1007ba337c44SJed Brown #undef __FUNCT__
1008ba337c44SJed Brown #define __FUNCT__ "MatConjugate_MPIDense"
10097087cfbeSBarry Smith PetscErrorCode  MatConjugate_MPIDense(Mat mat)
1010ba337c44SJed Brown {
1011ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense *)mat->data;
1012ba337c44SJed Brown   PetscErrorCode ierr;
1013ba337c44SJed Brown 
1014ba337c44SJed Brown   PetscFunctionBegin;
1015ba337c44SJed Brown   ierr = MatConjugate(a->A);CHKERRQ(ierr);
1016ba337c44SJed Brown   PetscFunctionReturn(0);
1017ba337c44SJed Brown }
1018ba337c44SJed Brown 
1019ba337c44SJed Brown #undef __FUNCT__
1020ba337c44SJed Brown #define __FUNCT__ "MatRealPart_MPIDense"
1021ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A)
1022ba337c44SJed Brown {
1023ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1024ba337c44SJed Brown   PetscErrorCode ierr;
1025ba337c44SJed Brown 
1026ba337c44SJed Brown   PetscFunctionBegin;
1027ba337c44SJed Brown   ierr = MatRealPart(a->A);CHKERRQ(ierr);
1028ba337c44SJed Brown   PetscFunctionReturn(0);
1029ba337c44SJed Brown }
1030ba337c44SJed Brown 
1031ba337c44SJed Brown #undef __FUNCT__
1032ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_MPIDense"
1033ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1034ba337c44SJed Brown {
1035ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1036ba337c44SJed Brown   PetscErrorCode ierr;
1037ba337c44SJed Brown 
1038ba337c44SJed Brown   PetscFunctionBegin;
1039ba337c44SJed Brown   ierr = MatImaginaryPart(a->A);CHKERRQ(ierr);
1040ba337c44SJed Brown   PetscFunctionReturn(0);
1041ba337c44SJed Brown }
1042ba337c44SJed Brown 
10430716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);
10440716a85fSBarry Smith #undef __FUNCT__
10450716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense"
10460716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
10470716a85fSBarry Smith {
10480716a85fSBarry Smith   PetscErrorCode ierr;
10490716a85fSBarry Smith   PetscInt       i,n;
10500716a85fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
10510716a85fSBarry Smith   PetscReal      *work;
10520716a85fSBarry Smith 
10530716a85fSBarry Smith   PetscFunctionBegin;
10540716a85fSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&n);CHKERRQ(ierr);
10550716a85fSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr);
10560716a85fSBarry Smith   ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr);
10570716a85fSBarry Smith   if (type == NORM_2) {
10580716a85fSBarry Smith     for (i=0; i<n; i++) work[i] *= work[i];
10590716a85fSBarry Smith   }
10600716a85fSBarry Smith   if (type == NORM_INFINITY) {
10610716a85fSBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr);
10620716a85fSBarry Smith   } else {
10630716a85fSBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr);
10640716a85fSBarry Smith   }
10650716a85fSBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
10660716a85fSBarry Smith   if (type == NORM_2) {
10678f1a2a5eSBarry Smith     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
10680716a85fSBarry Smith   }
10690716a85fSBarry Smith   PetscFunctionReturn(0);
10700716a85fSBarry Smith }
10710716a85fSBarry Smith 
107273a71a0fSBarry Smith #undef __FUNCT__
107373a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_MPIDense"
107473a71a0fSBarry Smith static PetscErrorCode  MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
107573a71a0fSBarry Smith {
107673a71a0fSBarry Smith   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;
107773a71a0fSBarry Smith   PetscErrorCode ierr;
107873a71a0fSBarry Smith   PetscScalar    *a;
107973a71a0fSBarry Smith   PetscInt       m,n,i;
108073a71a0fSBarry Smith 
108173a71a0fSBarry Smith   PetscFunctionBegin;
108273a71a0fSBarry Smith   ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr);
1083*8c778c55SBarry Smith   ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr);
108473a71a0fSBarry Smith   for (i=0; i<m*n; i++) {
108573a71a0fSBarry Smith     ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr);
108673a71a0fSBarry Smith   }
1087*8c778c55SBarry Smith   ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr);
108873a71a0fSBarry Smith   PetscFunctionReturn(0);
108973a71a0fSBarry Smith }
109073a71a0fSBarry Smith 
10918965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
109209dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
109309dc0095SBarry Smith        MatGetRow_MPIDense,
109409dc0095SBarry Smith        MatRestoreRow_MPIDense,
109509dc0095SBarry Smith        MatMult_MPIDense,
109697304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
10977c922b88SBarry Smith        MatMultTranspose_MPIDense,
10987c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
10998965ea79SLois Curfman McInnes        0,
110009dc0095SBarry Smith        0,
110109dc0095SBarry Smith        0,
110297304618SKris Buschelman /*10*/ 0,
110309dc0095SBarry Smith        0,
110409dc0095SBarry Smith        0,
110509dc0095SBarry Smith        0,
110609dc0095SBarry Smith        MatTranspose_MPIDense,
110797304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
11086e4ee0c6SHong Zhang        MatEqual_MPIDense,
110909dc0095SBarry Smith        MatGetDiagonal_MPIDense,
11105b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
111109dc0095SBarry Smith        MatNorm_MPIDense,
111297304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
111309dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
111409dc0095SBarry Smith        MatSetOption_MPIDense,
111509dc0095SBarry Smith        MatZeroEntries_MPIDense,
1116d519adbfSMatthew Knepley /*24*/ MatZeroRows_MPIDense,
1117919b68f7SBarry Smith        0,
111801b82886SBarry Smith        0,
111901b82886SBarry Smith        0,
112001b82886SBarry Smith        0,
11214994cf47SJed Brown /*29*/ MatSetUp_MPIDense,
1122273d9f13SBarry Smith        0,
112309dc0095SBarry Smith        0,
1124*8c778c55SBarry Smith        0,
1125*8c778c55SBarry Smith        0,
1126d519adbfSMatthew Knepley /*34*/ MatDuplicate_MPIDense,
112709dc0095SBarry Smith        0,
112809dc0095SBarry Smith        0,
112909dc0095SBarry Smith        0,
113009dc0095SBarry Smith        0,
1131d519adbfSMatthew Knepley /*39*/ MatAXPY_MPIDense,
11322ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
113309dc0095SBarry Smith        0,
113409dc0095SBarry Smith        MatGetValues_MPIDense,
113509dc0095SBarry Smith        0,
1136d519adbfSMatthew Knepley /*44*/ 0,
113709dc0095SBarry Smith        MatScale_MPIDense,
113809dc0095SBarry Smith        0,
113909dc0095SBarry Smith        0,
114009dc0095SBarry Smith        0,
114173a71a0fSBarry Smith /*49*/ MatSetRandom_MPIDense,
114209dc0095SBarry Smith        0,
114309dc0095SBarry Smith        0,
114409dc0095SBarry Smith        0,
114509dc0095SBarry Smith        0,
1146d519adbfSMatthew Knepley /*54*/ 0,
114709dc0095SBarry Smith        0,
114809dc0095SBarry Smith        0,
114909dc0095SBarry Smith        0,
115009dc0095SBarry Smith        0,
1151d519adbfSMatthew Knepley /*59*/ MatGetSubMatrix_MPIDense,
1152b9b97703SBarry Smith        MatDestroy_MPIDense,
1153b9b97703SBarry Smith        MatView_MPIDense,
1154357abbc8SBarry Smith        0,
115597304618SKris Buschelman        0,
1156d519adbfSMatthew Knepley /*64*/ 0,
115797304618SKris Buschelman        0,
115897304618SKris Buschelman        0,
115997304618SKris Buschelman        0,
116097304618SKris Buschelman        0,
1161d519adbfSMatthew Knepley /*69*/ 0,
116297304618SKris Buschelman        0,
116397304618SKris Buschelman        0,
116497304618SKris Buschelman        0,
116597304618SKris Buschelman        0,
1166d519adbfSMatthew Knepley /*74*/ 0,
116797304618SKris Buschelman        0,
116897304618SKris Buschelman        0,
116997304618SKris Buschelman        0,
117097304618SKris Buschelman        0,
1171d519adbfSMatthew Knepley /*79*/ 0,
117297304618SKris Buschelman        0,
117397304618SKris Buschelman        0,
117497304618SKris Buschelman        0,
11755bba2384SShri Abhyankar /*83*/ MatLoad_MPIDense,
1176865e5f61SKris Buschelman        0,
1177865e5f61SKris Buschelman        0,
1178865e5f61SKris Buschelman        0,
1179865e5f61SKris Buschelman        0,
1180865e5f61SKris Buschelman        0,
1181d519adbfSMatthew Knepley /*89*/
1182865e5f61SKris Buschelman        0,
1183865e5f61SKris Buschelman        0,
1184865e5f61SKris Buschelman        0,
11852fbe02b9SBarry Smith        0,
1186ba337c44SJed Brown        0,
1187d519adbfSMatthew Knepley /*94*/ 0,
1188865e5f61SKris Buschelman        0,
1189865e5f61SKris Buschelman        0,
1190ba337c44SJed Brown        0,
1191ba337c44SJed Brown        0,
1192ba337c44SJed Brown /*99*/ 0,
1193ba337c44SJed Brown        0,
1194ba337c44SJed Brown        0,
1195ba337c44SJed Brown        MatConjugate_MPIDense,
1196ba337c44SJed Brown        0,
1197ba337c44SJed Brown /*104*/0,
1198ba337c44SJed Brown        MatRealPart_MPIDense,
1199ba337c44SJed Brown        MatImaginaryPart_MPIDense,
120086d161a7SShri Abhyankar        0,
120186d161a7SShri Abhyankar        0,
120286d161a7SShri Abhyankar /*109*/0,
120386d161a7SShri Abhyankar        0,
120486d161a7SShri Abhyankar        0,
120586d161a7SShri Abhyankar        0,
120686d161a7SShri Abhyankar        0,
120786d161a7SShri Abhyankar /*114*/0,
120886d161a7SShri Abhyankar        0,
120986d161a7SShri Abhyankar        0,
121086d161a7SShri Abhyankar        0,
121186d161a7SShri Abhyankar        0,
121286d161a7SShri Abhyankar /*119*/0,
121386d161a7SShri Abhyankar        0,
121486d161a7SShri Abhyankar        0,
12150716a85fSBarry Smith        0,
12160716a85fSBarry Smith        0,
12170716a85fSBarry Smith /*124*/0,
12180716a85fSBarry Smith        MatGetColumnNorms_MPIDense
1219ba337c44SJed Brown };
12208965ea79SLois Curfman McInnes 
1221273d9f13SBarry Smith EXTERN_C_BEGIN
12224a2ae208SSatish Balay #undef __FUNCT__
1223a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
12247087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1225a23d5eceSKris Buschelman {
1226a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1227dfbe8321SBarry Smith   PetscErrorCode ierr;
1228a23d5eceSKris Buschelman 
1229a23d5eceSKris Buschelman   PetscFunctionBegin;
1230a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1231a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1232a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1233a23d5eceSKris Buschelman 
1234a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
123534ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
123634ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
123734ef9618SShri Abhyankar   a->nvec = mat->cmap->n;
123834ef9618SShri Abhyankar 
1239f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1240d0f46423SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr);
12415c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
12425c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
124352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1244a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1245a23d5eceSKris Buschelman }
1246a23d5eceSKris Buschelman EXTERN_C_END
1247a23d5eceSKris Buschelman 
1248a23d5eceSKris Buschelman EXTERN_C_BEGIN
1249a23d5eceSKris Buschelman #undef __FUNCT__
12504a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
12517087cfbeSBarry Smith PetscErrorCode  MatCreate_MPIDense(Mat mat)
1252273d9f13SBarry Smith {
1253273d9f13SBarry Smith   Mat_MPIDense   *a;
1254dfbe8321SBarry Smith   PetscErrorCode ierr;
1255273d9f13SBarry Smith 
1256273d9f13SBarry Smith   PetscFunctionBegin;
125738f2d2fdSLisandro Dalcin   ierr              = PetscNewLog(mat,Mat_MPIDense,&a);CHKERRQ(ierr);
1258b0a32e0cSBarry Smith   mat->data         = (void*)a;
1259273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1260273d9f13SBarry Smith 
1261273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
12627adad957SLisandro Dalcin   ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&a->rank);CHKERRQ(ierr);
12637adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)mat)->comm,&a->size);CHKERRQ(ierr);
1264273d9f13SBarry Smith 
1265273d9f13SBarry Smith   /* build cache for off array entries formed */
1266273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
12677adad957SLisandro Dalcin   ierr = MatStashCreate_Private(((PetscObject)mat)->comm,1,&mat->stash);CHKERRQ(ierr);
1268273d9f13SBarry Smith 
1269273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1270273d9f13SBarry Smith   a->lvec        = 0;
1271273d9f13SBarry Smith   a->Mvctx       = 0;
1272273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1273273d9f13SBarry Smith 
1274*8c778c55SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseGetArray_C","MatDenseGetArray_MPIDense",MatDenseGetArray_MPIDense);CHKERRQ(ierr);
1275*8c778c55SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseRestoreArray_C","MatDenseRestoreArray_MPIDense",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr);
1276*8c778c55SBarry Smith 
1277273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1278273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1279273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1280a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1281a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1282a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
12834ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",
12844ae313f4SHong Zhang                                      "MatMatMult_MPIAIJ_MPIDense",
12854ae313f4SHong Zhang                                       MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
12864ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",
12874ae313f4SHong Zhang                                      "MatMatMultSymbolic_MPIAIJ_MPIDense",
12884ae313f4SHong Zhang                                       MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
12894ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",
12904ae313f4SHong Zhang                                      "MatMatMultNumeric_MPIAIJ_MPIDense",
12914ae313f4SHong Zhang                                       MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
129238aed534SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr);
129301b82886SBarry Smith 
1294273d9f13SBarry Smith   PetscFunctionReturn(0);
1295273d9f13SBarry Smith }
1296273d9f13SBarry Smith EXTERN_C_END
1297273d9f13SBarry Smith 
1298209238afSKris Buschelman /*MC
1299002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1300209238afSKris Buschelman 
1301209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1302209238afSKris Buschelman    and MATMPIDENSE otherwise.
1303209238afSKris Buschelman 
1304209238afSKris Buschelman    Options Database Keys:
1305209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1306209238afSKris Buschelman 
1307209238afSKris Buschelman   Level: beginner
1308209238afSKris Buschelman 
130901b82886SBarry Smith 
1310209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1311209238afSKris Buschelman M*/
1312209238afSKris Buschelman 
13134a2ae208SSatish Balay #undef __FUNCT__
13144a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1315273d9f13SBarry Smith /*@C
1316273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1317273d9f13SBarry Smith 
1318273d9f13SBarry Smith    Not collective
1319273d9f13SBarry Smith 
1320273d9f13SBarry Smith    Input Parameters:
1321273d9f13SBarry Smith .  A - the matrix
1322273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1323273d9f13SBarry Smith    to control all matrix memory allocation.
1324273d9f13SBarry Smith 
1325273d9f13SBarry Smith    Notes:
1326273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1327273d9f13SBarry Smith    storage by columns.
1328273d9f13SBarry Smith 
1329273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1330273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1331273d9f13SBarry Smith    set data=PETSC_NULL.
1332273d9f13SBarry Smith 
1333273d9f13SBarry Smith    Level: intermediate
1334273d9f13SBarry Smith 
1335273d9f13SBarry Smith .keywords: matrix,dense, parallel
1336273d9f13SBarry Smith 
1337273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1338273d9f13SBarry Smith @*/
13397087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1340273d9f13SBarry Smith {
13414ac538c5SBarry Smith   PetscErrorCode ierr;
1342273d9f13SBarry Smith 
1343273d9f13SBarry Smith   PetscFunctionBegin;
13444ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar *),(mat,data));CHKERRQ(ierr);
1345273d9f13SBarry Smith   PetscFunctionReturn(0);
1346273d9f13SBarry Smith }
1347273d9f13SBarry Smith 
13484a2ae208SSatish Balay #undef __FUNCT__
134969b1f4b7SBarry Smith #define __FUNCT__ "MatCreateDense"
13508965ea79SLois Curfman McInnes /*@C
135169b1f4b7SBarry Smith    MatCreateDense - Creates a parallel matrix in dense format.
13528965ea79SLois Curfman McInnes 
1353db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1354db81eaa0SLois Curfman McInnes 
13558965ea79SLois Curfman McInnes    Input Parameters:
1356db81eaa0SLois Curfman McInnes +  comm - MPI communicator
13578965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1358db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
13598965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1360db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
13617f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1362dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
13638965ea79SLois Curfman McInnes 
13648965ea79SLois Curfman McInnes    Output Parameter:
1365477f1c0bSLois Curfman McInnes .  A - the matrix
13668965ea79SLois Curfman McInnes 
1367b259b22eSLois Curfman McInnes    Notes:
136839ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
136939ddd567SLois Curfman McInnes    storage by columns.
13708965ea79SLois Curfman McInnes 
137118f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
137218f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
13737f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
137418f449edSLois Curfman McInnes 
13758965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
13768965ea79SLois Curfman McInnes    (possibly both).
13778965ea79SLois Curfman McInnes 
1378027ccd11SLois Curfman McInnes    Level: intermediate
1379027ccd11SLois Curfman McInnes 
138039ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
13818965ea79SLois Curfman McInnes 
138239ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
13838965ea79SLois Curfman McInnes @*/
138469b1f4b7SBarry Smith PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
13858965ea79SLois Curfman McInnes {
13866849ba73SBarry Smith   PetscErrorCode ierr;
138713f74950SBarry Smith   PetscMPIInt    size;
13888965ea79SLois Curfman McInnes 
13893a40ed3dSBarry Smith   PetscFunctionBegin;
1390f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1391f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1392273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1393273d9f13SBarry Smith   if (size > 1) {
1394273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1395273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1396273d9f13SBarry Smith   } else {
1397273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1398273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
13998c469469SLois Curfman McInnes   }
14003a40ed3dSBarry Smith   PetscFunctionReturn(0);
14018965ea79SLois Curfman McInnes }
14028965ea79SLois Curfman McInnes 
14034a2ae208SSatish Balay #undef __FUNCT__
14044a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
14056849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
14068965ea79SLois Curfman McInnes {
14078965ea79SLois Curfman McInnes   Mat            mat;
14083501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1409dfbe8321SBarry Smith   PetscErrorCode ierr;
14108965ea79SLois Curfman McInnes 
14113a40ed3dSBarry Smith   PetscFunctionBegin;
14128965ea79SLois Curfman McInnes   *newmat       = 0;
14137adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&mat);CHKERRQ(ierr);
1414d0f46423SBarry Smith   ierr = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
14157adad957SLisandro Dalcin   ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1416834f8fabSBarry Smith   a                 = (Mat_MPIDense*)mat->data;
1417e04c1aa4SHong Zhang   ierr              = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
14185aa7edbeSHong Zhang 
1419d5f3da31SBarry Smith   mat->factortype   = A->factortype;
1420c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1421273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
14228965ea79SLois Curfman McInnes 
14238965ea79SLois Curfman McInnes   a->size           = oldmat->size;
14248965ea79SLois Curfman McInnes   a->rank           = oldmat->rank;
1425e0fa3b82SLois Curfman McInnes   mat->insertmode   = NOT_SET_VALUES;
1426b9b97703SBarry Smith   a->nvec           = oldmat->nvec;
14273782ba37SSatish Balay   a->donotstash     = oldmat->donotstash;
1428e04c1aa4SHong Zhang 
14291e1e43feSBarry Smith   ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr);
14301e1e43feSBarry Smith   ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr);
14318965ea79SLois Curfman McInnes 
1432329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
14335609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
143452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
143501b82886SBarry Smith 
14368965ea79SLois Curfman McInnes   *newmat = mat;
14373a40ed3dSBarry Smith   PetscFunctionReturn(0);
14388965ea79SLois Curfman McInnes }
14398965ea79SLois Curfman McInnes 
14404a2ae208SSatish Balay #undef __FUNCT__
14415bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
14425bba2384SShri Abhyankar PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat,PetscInt sizesset)
144386d161a7SShri Abhyankar {
144486d161a7SShri Abhyankar   PetscErrorCode ierr;
144586d161a7SShri Abhyankar   PetscMPIInt    rank,size;
144686d161a7SShri Abhyankar   PetscInt       *rowners,i,m,nz,j;
144786d161a7SShri Abhyankar   PetscScalar    *array,*vals,*vals_ptr;
144886d161a7SShri Abhyankar 
144986d161a7SShri Abhyankar   PetscFunctionBegin;
145086d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
145186d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
145286d161a7SShri Abhyankar 
145386d161a7SShri Abhyankar   /* determine ownership of all rows */
145486d161a7SShri Abhyankar   if (newmat->rmap->n < 0) m          = M/size + ((M % size) > rank);
145586d161a7SShri Abhyankar   else m = newmat->rmap->n;
145686d161a7SShri Abhyankar   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
145786d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
145886d161a7SShri Abhyankar   rowners[0] = 0;
145986d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
146086d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
146186d161a7SShri Abhyankar   }
146286d161a7SShri Abhyankar 
146386d161a7SShri Abhyankar   if (!sizesset) {
146486d161a7SShri Abhyankar     ierr = MatSetSizes(newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
146586d161a7SShri Abhyankar   }
146686d161a7SShri Abhyankar   ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1467*8c778c55SBarry Smith   ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr);
146886d161a7SShri Abhyankar 
146986d161a7SShri Abhyankar   if (!rank) {
147086d161a7SShri Abhyankar     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
147186d161a7SShri Abhyankar 
147286d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
147386d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
147486d161a7SShri Abhyankar 
147586d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
147686d161a7SShri Abhyankar     vals_ptr = vals;
147786d161a7SShri Abhyankar     for (i=0; i<m; i++) {
147886d161a7SShri Abhyankar       for (j=0; j<N; j++) {
147986d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
148086d161a7SShri Abhyankar       }
148186d161a7SShri Abhyankar     }
148286d161a7SShri Abhyankar 
148386d161a7SShri Abhyankar     /* read in other processors and ship out */
148486d161a7SShri Abhyankar     for (i=1; i<size; i++) {
148586d161a7SShri Abhyankar       nz   = (rowners[i+1] - rowners[i])*N;
148686d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1487a25532f0SBarry Smith       ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
148886d161a7SShri Abhyankar     }
148986d161a7SShri Abhyankar   } else {
149086d161a7SShri Abhyankar     /* receive numeric values */
149186d161a7SShri Abhyankar     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
149286d161a7SShri Abhyankar 
149386d161a7SShri Abhyankar     /* receive message of values*/
1494a25532f0SBarry Smith     ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
149586d161a7SShri Abhyankar 
149686d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
149786d161a7SShri Abhyankar     vals_ptr = vals;
149886d161a7SShri Abhyankar     for (i=0; i<m; i++) {
149986d161a7SShri Abhyankar       for (j=0; j<N; j++) {
150086d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
150186d161a7SShri Abhyankar       }
150286d161a7SShri Abhyankar     }
150386d161a7SShri Abhyankar   }
1504*8c778c55SBarry Smith   ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr);
150586d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
150686d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
150786d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150886d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150986d161a7SShri Abhyankar   PetscFunctionReturn(0);
151086d161a7SShri Abhyankar }
151186d161a7SShri Abhyankar 
151286d161a7SShri Abhyankar #undef __FUNCT__
15135bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense"
1514112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer)
151586d161a7SShri Abhyankar {
151686d161a7SShri Abhyankar   PetscScalar    *vals,*svals;
151786d161a7SShri Abhyankar   MPI_Comm       comm = ((PetscObject)viewer)->comm;
151886d161a7SShri Abhyankar   MPI_Status     status;
151986d161a7SShri Abhyankar   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
152086d161a7SShri Abhyankar   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
152186d161a7SShri Abhyankar   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
152286d161a7SShri Abhyankar   PetscInt       i,nz,j,rstart,rend,sizesset=1,grows,gcols;
152386d161a7SShri Abhyankar   int            fd;
152486d161a7SShri Abhyankar   PetscErrorCode ierr;
152586d161a7SShri Abhyankar 
152686d161a7SShri Abhyankar   PetscFunctionBegin;
152786d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
152886d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
152986d161a7SShri Abhyankar   if (!rank) {
153086d161a7SShri Abhyankar     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
153186d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
153286d161a7SShri Abhyankar     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
153386d161a7SShri Abhyankar   }
153486d161a7SShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) sizesset = 0;
153586d161a7SShri Abhyankar 
153686d161a7SShri Abhyankar   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
153786d161a7SShri Abhyankar   M = header[1]; N = header[2]; nz = header[3];
153886d161a7SShri Abhyankar 
153986d161a7SShri Abhyankar   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
154086d161a7SShri Abhyankar   if (sizesset && newmat->rmap->N < 0) newmat->rmap->N = M;
154186d161a7SShri Abhyankar   if (sizesset && newmat->cmap->N < 0) newmat->cmap->N = N;
154286d161a7SShri Abhyankar 
154386d161a7SShri Abhyankar   /* If global sizes are set, check if they are consistent with that given in the file */
154486d161a7SShri Abhyankar   if (sizesset) {
154586d161a7SShri Abhyankar     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
154686d161a7SShri Abhyankar   }
1547abd38a8fSBarry Smith   if (sizesset && newmat->rmap->N != grows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%d) and input matrix has (%d)",M,grows);
1548abd38a8fSBarry Smith   if (sizesset && newmat->cmap->N != gcols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%d) and input matrix has (%d)",N,gcols);
154986d161a7SShri Abhyankar 
155086d161a7SShri Abhyankar   /*
155186d161a7SShri Abhyankar        Handle case where matrix is stored on disk as a dense matrix
155286d161a7SShri Abhyankar   */
155386d161a7SShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
15545bba2384SShri Abhyankar     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat,sizesset);CHKERRQ(ierr);
155586d161a7SShri Abhyankar     PetscFunctionReturn(0);
155686d161a7SShri Abhyankar   }
155786d161a7SShri Abhyankar 
155886d161a7SShri Abhyankar   /* determine ownership of all rows */
155986d161a7SShri Abhyankar   if (newmat->rmap->n < 0) m          = PetscMPIIntCast(M/size + ((M % size) > rank));
156086d161a7SShri Abhyankar   else m = PetscMPIIntCast(newmat->rmap->n);
156186d161a7SShri Abhyankar   ierr       = PetscMalloc((size+2)*sizeof(PetscMPIInt),&rowners);CHKERRQ(ierr);
156286d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
156386d161a7SShri Abhyankar   rowners[0] = 0;
156486d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
156586d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
156686d161a7SShri Abhyankar   }
156786d161a7SShri Abhyankar   rstart = rowners[rank];
156886d161a7SShri Abhyankar   rend   = rowners[rank+1];
156986d161a7SShri Abhyankar 
157086d161a7SShri Abhyankar   /* distribute row lengths to all processors */
157186d161a7SShri Abhyankar   ierr    = PetscMalloc2(rend-rstart,PetscInt,&ourlens,rend-rstart,PetscInt,&offlens);CHKERRQ(ierr);
157286d161a7SShri Abhyankar   if (!rank) {
157386d161a7SShri Abhyankar     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
157486d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
157586d161a7SShri Abhyankar     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
157686d161a7SShri Abhyankar     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
157786d161a7SShri Abhyankar     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
157886d161a7SShri Abhyankar     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
157986d161a7SShri Abhyankar   } else {
158086d161a7SShri Abhyankar     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
158186d161a7SShri Abhyankar   }
158286d161a7SShri Abhyankar 
158386d161a7SShri Abhyankar   if (!rank) {
158486d161a7SShri Abhyankar     /* calculate the number of nonzeros on each processor */
158586d161a7SShri Abhyankar     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
158686d161a7SShri Abhyankar     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
158786d161a7SShri Abhyankar     for (i=0; i<size; i++) {
158886d161a7SShri Abhyankar       for (j=rowners[i]; j< rowners[i+1]; j++) {
158986d161a7SShri Abhyankar         procsnz[i] += rowlengths[j];
159086d161a7SShri Abhyankar       }
159186d161a7SShri Abhyankar     }
159286d161a7SShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
159386d161a7SShri Abhyankar 
159486d161a7SShri Abhyankar     /* determine max buffer needed and allocate it */
159586d161a7SShri Abhyankar     maxnz = 0;
159686d161a7SShri Abhyankar     for (i=0; i<size; i++) {
159786d161a7SShri Abhyankar       maxnz = PetscMax(maxnz,procsnz[i]);
159886d161a7SShri Abhyankar     }
159986d161a7SShri Abhyankar     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
160086d161a7SShri Abhyankar 
160186d161a7SShri Abhyankar     /* read in my part of the matrix column indices  */
160286d161a7SShri Abhyankar     nz = procsnz[0];
160386d161a7SShri Abhyankar     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
160486d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
160586d161a7SShri Abhyankar 
160686d161a7SShri Abhyankar     /* read in every one elses and ship off */
160786d161a7SShri Abhyankar     for (i=1; i<size; i++) {
160886d161a7SShri Abhyankar       nz   = procsnz[i];
160986d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
161086d161a7SShri Abhyankar       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
161186d161a7SShri Abhyankar     }
161286d161a7SShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
161386d161a7SShri Abhyankar   } else {
161486d161a7SShri Abhyankar     /* determine buffer space needed for message */
161586d161a7SShri Abhyankar     nz = 0;
161686d161a7SShri Abhyankar     for (i=0; i<m; i++) {
161786d161a7SShri Abhyankar       nz += ourlens[i];
161886d161a7SShri Abhyankar     }
161986d161a7SShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
162086d161a7SShri Abhyankar 
162186d161a7SShri Abhyankar     /* receive message of column indices*/
162286d161a7SShri Abhyankar     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
162386d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
162486d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
162586d161a7SShri Abhyankar   }
162686d161a7SShri Abhyankar 
162786d161a7SShri Abhyankar   /* loop over local rows, determining number of off diagonal entries */
162886d161a7SShri Abhyankar   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
162986d161a7SShri Abhyankar   jj = 0;
163086d161a7SShri Abhyankar   for (i=0; i<m; i++) {
163186d161a7SShri Abhyankar     for (j=0; j<ourlens[i]; j++) {
163286d161a7SShri Abhyankar       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
163386d161a7SShri Abhyankar       jj++;
163486d161a7SShri Abhyankar     }
163586d161a7SShri Abhyankar   }
163686d161a7SShri Abhyankar 
163786d161a7SShri Abhyankar   /* create our matrix */
163886d161a7SShri Abhyankar   for (i=0; i<m; i++) {
163986d161a7SShri Abhyankar     ourlens[i] -= offlens[i];
164086d161a7SShri Abhyankar   }
164186d161a7SShri Abhyankar 
164286d161a7SShri Abhyankar   if (!sizesset) {
164386d161a7SShri Abhyankar     ierr = MatSetSizes(newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
164486d161a7SShri Abhyankar   }
164586d161a7SShri Abhyankar   ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
164686d161a7SShri Abhyankar   for (i=0; i<m; i++) {
164786d161a7SShri Abhyankar     ourlens[i] += offlens[i];
164886d161a7SShri Abhyankar   }
164986d161a7SShri Abhyankar 
165086d161a7SShri Abhyankar   if (!rank) {
165186d161a7SShri Abhyankar     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
165286d161a7SShri Abhyankar 
165386d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
165486d161a7SShri Abhyankar     nz = procsnz[0];
165586d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
165686d161a7SShri Abhyankar 
165786d161a7SShri Abhyankar     /* insert into matrix */
165886d161a7SShri Abhyankar     jj      = rstart;
165986d161a7SShri Abhyankar     smycols = mycols;
166086d161a7SShri Abhyankar     svals   = vals;
166186d161a7SShri Abhyankar     for (i=0; i<m; i++) {
166286d161a7SShri Abhyankar       ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
166386d161a7SShri Abhyankar       smycols += ourlens[i];
166486d161a7SShri Abhyankar       svals   += ourlens[i];
166586d161a7SShri Abhyankar       jj++;
166686d161a7SShri Abhyankar     }
166786d161a7SShri Abhyankar 
166886d161a7SShri Abhyankar     /* read in other processors and ship out */
166986d161a7SShri Abhyankar     for (i=1; i<size; i++) {
167086d161a7SShri Abhyankar       nz   = procsnz[i];
167186d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
167286d161a7SShri Abhyankar       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr);
167386d161a7SShri Abhyankar     }
167486d161a7SShri Abhyankar     ierr = PetscFree(procsnz);CHKERRQ(ierr);
167586d161a7SShri Abhyankar   } else {
167686d161a7SShri Abhyankar     /* receive numeric values */
167786d161a7SShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
167886d161a7SShri Abhyankar 
167986d161a7SShri Abhyankar     /* receive message of values*/
168086d161a7SShri Abhyankar     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr);
168186d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
168286d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
168386d161a7SShri Abhyankar 
168486d161a7SShri Abhyankar     /* insert into matrix */
168586d161a7SShri Abhyankar     jj      = rstart;
168686d161a7SShri Abhyankar     smycols = mycols;
168786d161a7SShri Abhyankar     svals   = vals;
168886d161a7SShri Abhyankar     for (i=0; i<m; i++) {
168986d161a7SShri Abhyankar       ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
169086d161a7SShri Abhyankar       smycols += ourlens[i];
169186d161a7SShri Abhyankar       svals   += ourlens[i];
169286d161a7SShri Abhyankar       jj++;
169386d161a7SShri Abhyankar     }
169486d161a7SShri Abhyankar   }
169586d161a7SShri Abhyankar   ierr = PetscFree2(ourlens,offlens);CHKERRQ(ierr);
169686d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
169786d161a7SShri Abhyankar   ierr = PetscFree(mycols);CHKERRQ(ierr);
169886d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
169986d161a7SShri Abhyankar 
170086d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
170186d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
170286d161a7SShri Abhyankar   PetscFunctionReturn(0);
170386d161a7SShri Abhyankar }
170486d161a7SShri Abhyankar 
170586d161a7SShri Abhyankar #undef __FUNCT__
17066e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense"
1707ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
17086e4ee0c6SHong Zhang {
17096e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
17106e4ee0c6SHong Zhang   Mat            a,b;
1711ace3abfcSBarry Smith   PetscBool      flg;
17126e4ee0c6SHong Zhang   PetscErrorCode ierr;
171390ace30eSBarry Smith 
17146e4ee0c6SHong Zhang   PetscFunctionBegin;
17156e4ee0c6SHong Zhang   a = matA->A;
17166e4ee0c6SHong Zhang   b = matB->A;
17176e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
17187adad957SLisandro Dalcin   ierr = MPI_Allreduce(&flg,flag,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr);
17196e4ee0c6SHong Zhang   PetscFunctionReturn(0);
17206e4ee0c6SHong Zhang }
172190ace30eSBarry Smith 
1722