xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 11aeaf0a669438c21ebbbd10be78fc394d59507d)
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);
898c778c55SBarry Smith     ierr = MatDenseGetArray(mdn->A,&array);CHKERRQ(ierr);
9011bd1e4dSLisandro Dalcin     ierr = MatSeqDenseSetPreallocation(B,array+m*rstart);CHKERRQ(ierr);
918c778c55SBarry 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__
1688c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_MPIDense"
1698c778c55SBarry 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;
1758c778c55SBarry 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__
2448c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_MPIDense"
2458c778c55SBarry 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;
2518c778c55SBarry 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 */
267c3aae356SJed Brown   ierr = MPI_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,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;
561aa482453SBarry Smith #if defined(PETSC_USE_LOG)
562d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N);
5638965ea79SLois Curfman McInnes #endif
5648798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
5656bf464f9SBarry Smith   ierr = MatDestroy(&mdn->A);CHKERRQ(ierr);
5666bf464f9SBarry Smith   ierr = VecDestroy(&mdn->lvec);CHKERRQ(ierr);
5676bf464f9SBarry Smith   ierr = VecScatterDestroy(&mdn->Mvctx);CHKERRQ(ierr);
56801b82886SBarry Smith 
569bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
570dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
571901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
572901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5734ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5744ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5754ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5763a40ed3dSBarry Smith   PetscFunctionReturn(0);
5778965ea79SLois Curfman McInnes }
57839ddd567SLois Curfman McInnes 
5794a2ae208SSatish Balay #undef __FUNCT__
5804a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5816849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5828965ea79SLois Curfman McInnes {
58339ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
584dfbe8321SBarry Smith   PetscErrorCode    ierr;
585aa05aa95SBarry Smith   PetscViewerFormat format;
586aa05aa95SBarry Smith   int               fd;
587d0f46423SBarry Smith   PetscInt          header[4],mmax,N = mat->cmap->N,i,j,m,k;
588aa05aa95SBarry Smith   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
589578230a0SSatish Balay   PetscScalar       *work,*v,*vv;
590aa05aa95SBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;
5917056b6fcSBarry Smith 
5923a40ed3dSBarry Smith   PetscFunctionBegin;
59339ddd567SLois Curfman McInnes   if (mdn->size == 1) {
59439ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
595aa05aa95SBarry Smith   } else {
596aa05aa95SBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
5977adad957SLisandro Dalcin     ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&rank);CHKERRQ(ierr);
5987adad957SLisandro Dalcin     ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr);
599aa05aa95SBarry Smith 
600aa05aa95SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
601f4403165SShri Abhyankar     if (format == PETSC_VIEWER_NATIVE) {
602aa05aa95SBarry Smith 
603aa05aa95SBarry Smith       if (!rank) {
604aa05aa95SBarry Smith         /* store the matrix as a dense matrix */
6050700a824SBarry Smith         header[0] = MAT_FILE_CLASSID;
606d0f46423SBarry Smith         header[1] = mat->rmap->N;
607aa05aa95SBarry Smith         header[2] = N;
608aa05aa95SBarry Smith         header[3] = MATRIX_BINARY_FORMAT_DENSE;
609aa05aa95SBarry Smith         ierr = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
610aa05aa95SBarry Smith 
611aa05aa95SBarry Smith         /* get largest work array needed for transposing array */
612d0f46423SBarry Smith         mmax = mat->rmap->n;
613aa05aa95SBarry Smith         for (i=1; i<size; i++) {
614d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
6158965ea79SLois Curfman McInnes         }
616aa05aa95SBarry Smith         ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&work);CHKERRQ(ierr);
617aa05aa95SBarry Smith 
618aa05aa95SBarry Smith         /* write out local array, by rows */
619d0f46423SBarry Smith         m    = mat->rmap->n;
620aa05aa95SBarry Smith         v    = a->v;
621aa05aa95SBarry Smith         for (j=0; j<N; j++) {
622aa05aa95SBarry Smith           for (i=0; i<m; i++) {
623578230a0SSatish Balay             work[j + i*N] = *v++;
624aa05aa95SBarry Smith           }
625aa05aa95SBarry Smith         }
626aa05aa95SBarry Smith         ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
627aa05aa95SBarry Smith         /* get largest work array to receive messages from other processes, excludes process zero */
628aa05aa95SBarry Smith         mmax = 0;
629aa05aa95SBarry Smith         for (i=1; i<size; i++) {
630d0f46423SBarry Smith           mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]);
631aa05aa95SBarry Smith         }
632578230a0SSatish Balay         ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&vv);CHKERRQ(ierr);
633aa05aa95SBarry Smith         for (k = 1; k < size; k++) {
634f8009846SMatthew Knepley           v    = vv;
635d0f46423SBarry Smith           m    = mat->rmap->range[k+1] - mat->rmap->range[k];
636a25532f0SBarry Smith           ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,((PetscObject)mat)->comm);CHKERRQ(ierr);
637aa05aa95SBarry Smith 
638aa05aa95SBarry Smith           for (j = 0; j < N; j++) {
639aa05aa95SBarry Smith             for (i = 0; i < m; i++) {
640578230a0SSatish Balay               work[j + i*N] = *v++;
641aa05aa95SBarry Smith             }
642aa05aa95SBarry Smith           }
643aa05aa95SBarry Smith           ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
644aa05aa95SBarry Smith         }
645aa05aa95SBarry Smith         ierr = PetscFree(work);CHKERRQ(ierr);
646578230a0SSatish Balay         ierr = PetscFree(vv);CHKERRQ(ierr);
647aa05aa95SBarry Smith       } else {
648a25532f0SBarry Smith         ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,((PetscObject)mat)->comm);CHKERRQ(ierr);
649aa05aa95SBarry Smith       }
650eb3f98d2SBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE)");
651aa05aa95SBarry Smith   }
6523a40ed3dSBarry Smith   PetscFunctionReturn(0);
6538965ea79SLois Curfman McInnes }
6548965ea79SLois Curfman McInnes 
6554a2ae208SSatish Balay #undef __FUNCT__
6564a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
6576849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6588965ea79SLois Curfman McInnes {
65939ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
660dfbe8321SBarry Smith   PetscErrorCode    ierr;
66132dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
66219fd82e9SBarry Smith   PetscViewerType   vtype;
663ace3abfcSBarry Smith   PetscBool         iascii,isdraw;
664b0a32e0cSBarry Smith   PetscViewer       sviewer;
665f3ef73ceSBarry Smith   PetscViewerFormat format;
6668965ea79SLois Curfman McInnes 
6673a40ed3dSBarry Smith   PetscFunctionBegin;
668251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
669251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
67032077d6dSBarry Smith   if (iascii) {
671b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
672b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
673456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6744e220ebcSLois Curfman McInnes       MatInfo info;
675888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
6767b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
6777b23a99aSBarry 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);
678b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6797b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
6805d00a290SHong Zhang       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
6813a40ed3dSBarry Smith       PetscFunctionReturn(0);
682fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6833a40ed3dSBarry Smith       PetscFunctionReturn(0);
6848965ea79SLois Curfman McInnes     }
685f1af5d2fSBarry Smith   } else if (isdraw) {
686b0a32e0cSBarry Smith     PetscDraw  draw;
687ace3abfcSBarry Smith     PetscBool  isnull;
688f1af5d2fSBarry Smith 
689b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
690b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
691f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
692f1af5d2fSBarry Smith   }
69377ed5343SBarry Smith 
6948965ea79SLois Curfman McInnes   if (size == 1) {
69539ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6963a40ed3dSBarry Smith   } else {
6978965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
6988965ea79SLois Curfman McInnes     Mat         A;
699d0f46423SBarry Smith     PetscInt    M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz;
700ba8c8a56SBarry Smith     PetscInt    *cols;
701ba8c8a56SBarry Smith     PetscScalar *vals;
7028965ea79SLois Curfman McInnes 
7037adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)mat)->comm,&A);CHKERRQ(ierr);
7048965ea79SLois Curfman McInnes     if (!rank) {
705f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
7063a40ed3dSBarry Smith     } else {
707f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
7088965ea79SLois Curfman McInnes     }
7097adad957SLisandro Dalcin     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
710878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
711a2ea699eSBarry Smith     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);CHKERRQ(ierr);
71252e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
7138965ea79SLois Curfman McInnes 
71439ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
71539ddd567SLois Curfman McInnes        but it's quick for now */
71651022da4SBarry Smith     A->insertmode = INSERT_VALUES;
717d0f46423SBarry Smith     row = mat->rmap->rstart; m = mdn->A->rmap->n;
7188965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
719ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
720ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
721ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
72239ddd567SLois Curfman McInnes       row++;
7238965ea79SLois Curfman McInnes     }
7248965ea79SLois Curfman McInnes 
7256d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7266d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
727b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
728b9b97703SBarry Smith     if (!rank) {
7297566de4bSShri Abhyankar       ierr = PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr);
7307566de4bSShri Abhyankar       /* Set the type name to MATMPIDense so that the correct type can be printed out by PetscObjectPrintClassNamePrefixType() in MatView_SeqDense_ASCII()*/
7317566de4bSShri Abhyankar       PetscStrcpy(((PetscObject)((Mat_MPIDense*)(A->data))->A)->type_name,MATMPIDENSE);
7326831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
7338965ea79SLois Curfman McInnes     }
734b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
735b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7366bf464f9SBarry Smith     ierr = MatDestroy(&A);CHKERRQ(ierr);
7378965ea79SLois Curfman McInnes   }
7383a40ed3dSBarry Smith   PetscFunctionReturn(0);
7398965ea79SLois Curfman McInnes }
7408965ea79SLois Curfman McInnes 
7414a2ae208SSatish Balay #undef __FUNCT__
7424a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
743dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
7448965ea79SLois Curfman McInnes {
745dfbe8321SBarry Smith   PetscErrorCode ierr;
746ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw,issocket;
7478965ea79SLois Curfman McInnes 
748433994e6SBarry Smith   PetscFunctionBegin;
749251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
750251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
751251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr);
752251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
7530f5bd95cSBarry Smith 
75432077d6dSBarry Smith   if (iascii || issocket || isdraw) {
755f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7560f5bd95cSBarry Smith   } else if (isbinary) {
7573a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
758*11aeaf0aSBarry Smith   }
7593a40ed3dSBarry Smith   PetscFunctionReturn(0);
7608965ea79SLois Curfman McInnes }
7618965ea79SLois Curfman McInnes 
7624a2ae208SSatish Balay #undef __FUNCT__
7634a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
764dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7658965ea79SLois Curfman McInnes {
7663501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7673501a2bdSLois Curfman McInnes   Mat            mdn = mat->A;
768dfbe8321SBarry Smith   PetscErrorCode ierr;
769329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7708965ea79SLois Curfman McInnes 
7713a40ed3dSBarry Smith   PetscFunctionBegin;
7724e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
7734e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7744e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7754e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7768965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7774e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7784e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7794e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7804e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7814e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7828965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
783d9822059SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,((PetscObject)A)->comm);CHKERRQ(ierr);
7844e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7854e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7864e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7874e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7884e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7898965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
790d9822059SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
7914e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7924e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7934e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7944e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7954e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7968965ea79SLois Curfman McInnes   }
7974e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
7984e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
7994e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
8003a40ed3dSBarry Smith   PetscFunctionReturn(0);
8018965ea79SLois Curfman McInnes }
8028965ea79SLois Curfman McInnes 
8034a2ae208SSatish Balay #undef __FUNCT__
8044a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
805ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool  flg)
8068965ea79SLois Curfman McInnes {
80739ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
808dfbe8321SBarry Smith   PetscErrorCode ierr;
8098965ea79SLois Curfman McInnes 
8103a40ed3dSBarry Smith   PetscFunctionBegin;
81112c028f9SKris Buschelman   switch (op) {
812512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
81312c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
81412c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
8154e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
81612c028f9SKris Buschelman     break;
81712c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
8184e0d8c25SBarry Smith     a->roworiented = flg;
8194e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82012c028f9SKris Buschelman     break;
8214e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
82213fa8e87SLisandro Dalcin   case MAT_KEEP_NONZERO_PATTERN:
82312c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
824290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
82512c028f9SKris Buschelman     break;
82612c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
8274e0d8c25SBarry Smith     a->donotstash = flg;
82812c028f9SKris Buschelman     break;
82977e54ba9SKris Buschelman   case MAT_SYMMETRIC:
83077e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8319a4540c5SBarry Smith   case MAT_HERMITIAN:
8329a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
833600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
834290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
83577e54ba9SKris Buschelman     break;
83612c028f9SKris Buschelman   default:
837e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
8383a40ed3dSBarry Smith   }
8393a40ed3dSBarry Smith   PetscFunctionReturn(0);
8408965ea79SLois Curfman McInnes }
8418965ea79SLois Curfman McInnes 
8428965ea79SLois Curfman McInnes 
8434a2ae208SSatish Balay #undef __FUNCT__
8444a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
845dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8465b2fa520SLois Curfman McInnes {
8475b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8485b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
84987828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
850dfbe8321SBarry Smith   PetscErrorCode ierr;
851d0f46423SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n;
8525b2fa520SLois Curfman McInnes 
8535b2fa520SLois Curfman McInnes   PetscFunctionBegin;
85472d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8555b2fa520SLois Curfman McInnes   if (ll) {
85672d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
857e32f2f54SBarry Smith     if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8581ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8595b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8605b2fa520SLois Curfman McInnes       x = l[i];
8615b2fa520SLois Curfman McInnes       v = mat->v + i;
8625b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8635b2fa520SLois Curfman McInnes     }
8641ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
865efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8665b2fa520SLois Curfman McInnes   }
8675b2fa520SLois Curfman McInnes   if (rr) {
868175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
869e32f2f54SBarry Smith     if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
870ca9f406cSSatish Balay     ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
871ca9f406cSSatish Balay     ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8721ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8735b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8745b2fa520SLois Curfman McInnes       x = r[i];
8755b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8765b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
8775b2fa520SLois Curfman McInnes     }
8781ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
879efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8805b2fa520SLois Curfman McInnes   }
8815b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8825b2fa520SLois Curfman McInnes }
8835b2fa520SLois Curfman McInnes 
8844a2ae208SSatish Balay #undef __FUNCT__
8854a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
886dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
887096963f5SLois Curfman McInnes {
8883501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8893501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
890dfbe8321SBarry Smith   PetscErrorCode ierr;
89113f74950SBarry Smith   PetscInt       i,j;
892329f5518SBarry Smith   PetscReal      sum = 0.0;
89387828ca2SBarry Smith   PetscScalar    *v = mat->v;
8943501a2bdSLois Curfman McInnes 
8953a40ed3dSBarry Smith   PetscFunctionBegin;
8963501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
897064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
8983501a2bdSLois Curfman McInnes   } else {
8993501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
900d0f46423SBarry Smith       for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) {
901329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
9023501a2bdSLois Curfman McInnes       }
903d9822059SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
9048f1a2a5eSBarry Smith       *nrm = PetscSqrtReal(*nrm);
905dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr);
9063a40ed3dSBarry Smith     } else if (type == NORM_1) {
907329f5518SBarry Smith       PetscReal *tmp,*tmp2;
90874ed9c26SBarry Smith       ierr = PetscMalloc2(A->cmap->N,PetscReal,&tmp,A->cmap->N,PetscReal,&tmp2);CHKERRQ(ierr);
90974ed9c26SBarry Smith       ierr = PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
91074ed9c26SBarry Smith       ierr = PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
911064f8208SBarry Smith       *nrm = 0.0;
9123501a2bdSLois Curfman McInnes       v = mat->v;
913d0f46423SBarry Smith       for (j=0; j<mdn->A->cmap->n; j++) {
914d0f46423SBarry Smith         for (i=0; i<mdn->A->rmap->n; i++) {
91567e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
9163501a2bdSLois Curfman McInnes         }
9173501a2bdSLois Curfman McInnes       }
918d9822059SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
919d0f46423SBarry Smith       for (j=0; j<A->cmap->N; j++) {
920064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
9213501a2bdSLois Curfman McInnes       }
92274ed9c26SBarry Smith       ierr = PetscFree2(tmp,tmp);CHKERRQ(ierr);
923d0f46423SBarry Smith       ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr);
9243a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
925329f5518SBarry Smith       PetscReal ntemp;
9263501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
927d9822059SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,((PetscObject)A)->comm);CHKERRQ(ierr);
928e7e72b3dSBarry Smith     } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support for two norm");
9293501a2bdSLois Curfman McInnes   }
9303a40ed3dSBarry Smith   PetscFunctionReturn(0);
9313501a2bdSLois Curfman McInnes }
9323501a2bdSLois Curfman McInnes 
9334a2ae208SSatish Balay #undef __FUNCT__
9344a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
935fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
9363501a2bdSLois Curfman McInnes {
9373501a2bdSLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
9383501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
9393501a2bdSLois Curfman McInnes   Mat            B;
940d0f46423SBarry Smith   PetscInt       M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart;
9416849ba73SBarry Smith   PetscErrorCode ierr;
94213f74950SBarry Smith   PetscInt       j,i;
94387828ca2SBarry Smith   PetscScalar    *v;
9443501a2bdSLois Curfman McInnes 
9453a40ed3dSBarry Smith   PetscFunctionBegin;
946e7e72b3dSBarry Smith   if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Supports square matrix only in-place");
947fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
9487adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
949d0f46423SBarry Smith     ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr);
9507adad957SLisandro Dalcin     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
951878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
952fc4dec0aSBarry Smith   } else {
953fc4dec0aSBarry Smith     B = *matout;
954fc4dec0aSBarry Smith   }
9553501a2bdSLois Curfman McInnes 
956d0f46423SBarry Smith   m = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v;
9571acff37aSSatish Balay   ierr = PetscMalloc(m*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
9583501a2bdSLois Curfman McInnes   for (i=0; i<m; i++) rwork[i] = rstart + i;
9591acff37aSSatish Balay   for (j=0; j<n; j++) {
9603501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9613501a2bdSLois Curfman McInnes     v   += m;
9623501a2bdSLois Curfman McInnes   }
963606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9646d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9656d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
966815cbec1SBarry Smith   if (reuse == MAT_INITIAL_MATRIX || *matout != A) {
9673501a2bdSLois Curfman McInnes     *matout = B;
9683501a2bdSLois Curfman McInnes   } else {
969eb6b5d47SBarry Smith     ierr = MatHeaderMerge(A,B);CHKERRQ(ierr);
9703501a2bdSLois Curfman McInnes   }
9713a40ed3dSBarry Smith   PetscFunctionReturn(0);
972096963f5SLois Curfman McInnes }
973096963f5SLois Curfman McInnes 
97444cd7ae7SLois Curfman McInnes 
9756849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
976d84338a6SBarry Smith extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar);
9778965ea79SLois Curfman McInnes 
9784a2ae208SSatish Balay #undef __FUNCT__
9794994cf47SJed Brown #define __FUNCT__ "MatSetUp_MPIDense"
9804994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A)
981273d9f13SBarry Smith {
982dfbe8321SBarry Smith   PetscErrorCode ierr;
983273d9f13SBarry Smith 
984273d9f13SBarry Smith   PetscFunctionBegin;
985273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
986273d9f13SBarry Smith   PetscFunctionReturn(0);
987273d9f13SBarry Smith }
988273d9f13SBarry Smith 
98930716080SHong Zhang #undef __FUNCT__
990488007eeSBarry Smith #define __FUNCT__ "MatAXPY_MPIDense"
991488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
992488007eeSBarry Smith {
993488007eeSBarry Smith   PetscErrorCode ierr;
994488007eeSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data;
995488007eeSBarry Smith 
996488007eeSBarry Smith   PetscFunctionBegin;
997488007eeSBarry Smith   ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr);
998488007eeSBarry Smith   PetscFunctionReturn(0);
999488007eeSBarry Smith }
1000488007eeSBarry Smith 
1001ba337c44SJed Brown #undef __FUNCT__
1002ba337c44SJed Brown #define __FUNCT__ "MatConjugate_MPIDense"
10037087cfbeSBarry Smith PetscErrorCode  MatConjugate_MPIDense(Mat mat)
1004ba337c44SJed Brown {
1005ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense *)mat->data;
1006ba337c44SJed Brown   PetscErrorCode ierr;
1007ba337c44SJed Brown 
1008ba337c44SJed Brown   PetscFunctionBegin;
1009ba337c44SJed Brown   ierr = MatConjugate(a->A);CHKERRQ(ierr);
1010ba337c44SJed Brown   PetscFunctionReturn(0);
1011ba337c44SJed Brown }
1012ba337c44SJed Brown 
1013ba337c44SJed Brown #undef __FUNCT__
1014ba337c44SJed Brown #define __FUNCT__ "MatRealPart_MPIDense"
1015ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A)
1016ba337c44SJed Brown {
1017ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1018ba337c44SJed Brown   PetscErrorCode ierr;
1019ba337c44SJed Brown 
1020ba337c44SJed Brown   PetscFunctionBegin;
1021ba337c44SJed Brown   ierr = MatRealPart(a->A);CHKERRQ(ierr);
1022ba337c44SJed Brown   PetscFunctionReturn(0);
1023ba337c44SJed Brown }
1024ba337c44SJed Brown 
1025ba337c44SJed Brown #undef __FUNCT__
1026ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_MPIDense"
1027ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1028ba337c44SJed Brown {
1029ba337c44SJed Brown   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
1030ba337c44SJed Brown   PetscErrorCode ierr;
1031ba337c44SJed Brown 
1032ba337c44SJed Brown   PetscFunctionBegin;
1033ba337c44SJed Brown   ierr = MatImaginaryPart(a->A);CHKERRQ(ierr);
1034ba337c44SJed Brown   PetscFunctionReturn(0);
1035ba337c44SJed Brown }
1036ba337c44SJed Brown 
10370716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*);
10380716a85fSBarry Smith #undef __FUNCT__
10390716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense"
10400716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
10410716a85fSBarry Smith {
10420716a85fSBarry Smith   PetscErrorCode ierr;
10430716a85fSBarry Smith   PetscInt       i,n;
10440716a85fSBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
10450716a85fSBarry Smith   PetscReal      *work;
10460716a85fSBarry Smith 
10470716a85fSBarry Smith   PetscFunctionBegin;
10480716a85fSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&n);CHKERRQ(ierr);
10490716a85fSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr);
10500716a85fSBarry Smith   ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr);
10510716a85fSBarry Smith   if (type == NORM_2) {
10520716a85fSBarry Smith     for (i=0; i<n; i++) work[i] *= work[i];
10530716a85fSBarry Smith   }
10540716a85fSBarry Smith   if (type == NORM_INFINITY) {
10550716a85fSBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr);
10560716a85fSBarry Smith   } else {
10570716a85fSBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr);
10580716a85fSBarry Smith   }
10590716a85fSBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
10600716a85fSBarry Smith   if (type == NORM_2) {
10618f1a2a5eSBarry Smith     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
10620716a85fSBarry Smith   }
10630716a85fSBarry Smith   PetscFunctionReturn(0);
10640716a85fSBarry Smith }
10650716a85fSBarry Smith 
106673a71a0fSBarry Smith #undef __FUNCT__
106773a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_MPIDense"
106873a71a0fSBarry Smith static PetscErrorCode  MatSetRandom_MPIDense(Mat x,PetscRandom rctx)
106973a71a0fSBarry Smith {
107073a71a0fSBarry Smith   Mat_MPIDense   *d = (Mat_MPIDense*)x->data;
107173a71a0fSBarry Smith   PetscErrorCode ierr;
107273a71a0fSBarry Smith   PetscScalar    *a;
107373a71a0fSBarry Smith   PetscInt       m,n,i;
107473a71a0fSBarry Smith 
107573a71a0fSBarry Smith   PetscFunctionBegin;
107673a71a0fSBarry Smith   ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr);
10778c778c55SBarry Smith   ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr);
107873a71a0fSBarry Smith   for (i=0; i<m*n; i++) {
107973a71a0fSBarry Smith     ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr);
108073a71a0fSBarry Smith   }
10818c778c55SBarry Smith   ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr);
108273a71a0fSBarry Smith   PetscFunctionReturn(0);
108373a71a0fSBarry Smith }
108473a71a0fSBarry Smith 
10858965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
108609dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense,
108709dc0095SBarry Smith                                         MatGetRow_MPIDense,
108809dc0095SBarry Smith                                         MatRestoreRow_MPIDense,
108909dc0095SBarry Smith                                         MatMult_MPIDense,
109097304618SKris Buschelman                                 /*  4*/ MatMultAdd_MPIDense,
10917c922b88SBarry Smith                                         MatMultTranspose_MPIDense,
10927c922b88SBarry Smith                                         MatMultTransposeAdd_MPIDense,
10938965ea79SLois Curfman McInnes                                         0,
109409dc0095SBarry Smith                                         0,
109509dc0095SBarry Smith                                         0,
109697304618SKris Buschelman                                 /* 10*/ 0,
109709dc0095SBarry Smith                                         0,
109809dc0095SBarry Smith                                         0,
109909dc0095SBarry Smith                                         0,
110009dc0095SBarry Smith                                         MatTranspose_MPIDense,
110197304618SKris Buschelman                                 /* 15*/ MatGetInfo_MPIDense,
11026e4ee0c6SHong Zhang                                         MatEqual_MPIDense,
110309dc0095SBarry Smith                                         MatGetDiagonal_MPIDense,
11045b2fa520SLois Curfman McInnes                                         MatDiagonalScale_MPIDense,
110509dc0095SBarry Smith                                         MatNorm_MPIDense,
110697304618SKris Buschelman                                 /* 20*/ MatAssemblyBegin_MPIDense,
110709dc0095SBarry Smith                                         MatAssemblyEnd_MPIDense,
110809dc0095SBarry Smith                                         MatSetOption_MPIDense,
110909dc0095SBarry Smith                                         MatZeroEntries_MPIDense,
1110d519adbfSMatthew Knepley                                 /* 24*/ MatZeroRows_MPIDense,
1111919b68f7SBarry Smith                                         0,
111201b82886SBarry Smith                                         0,
111301b82886SBarry Smith                                         0,
111401b82886SBarry Smith                                         0,
11154994cf47SJed Brown                                 /* 29*/ MatSetUp_MPIDense,
1116273d9f13SBarry Smith                                         0,
111709dc0095SBarry Smith                                         0,
11188c778c55SBarry Smith                                         0,
11198c778c55SBarry Smith                                         0,
1120d519adbfSMatthew Knepley                                 /* 34*/ MatDuplicate_MPIDense,
112109dc0095SBarry Smith                                         0,
112209dc0095SBarry Smith                                         0,
112309dc0095SBarry Smith                                         0,
112409dc0095SBarry Smith                                         0,
1125d519adbfSMatthew Knepley                                 /* 39*/ MatAXPY_MPIDense,
11262ce60cd0SSatish Balay                                         MatGetSubMatrices_MPIDense,
112709dc0095SBarry Smith                                         0,
112809dc0095SBarry Smith                                         MatGetValues_MPIDense,
112909dc0095SBarry Smith                                         0,
1130d519adbfSMatthew Knepley                                 /* 44*/ 0,
113109dc0095SBarry Smith                                         MatScale_MPIDense,
113209dc0095SBarry Smith                                         0,
113309dc0095SBarry Smith                                         0,
113409dc0095SBarry Smith                                         0,
113573a71a0fSBarry Smith                                 /* 49*/ MatSetRandom_MPIDense,
113609dc0095SBarry Smith                                         0,
113709dc0095SBarry Smith                                         0,
113809dc0095SBarry Smith                                         0,
113909dc0095SBarry Smith                                         0,
1140d519adbfSMatthew Knepley                                 /* 54*/ 0,
114109dc0095SBarry Smith                                         0,
114209dc0095SBarry Smith                                         0,
114309dc0095SBarry Smith                                         0,
114409dc0095SBarry Smith                                         0,
1145d519adbfSMatthew Knepley                                 /* 59*/ MatGetSubMatrix_MPIDense,
1146b9b97703SBarry Smith                                         MatDestroy_MPIDense,
1147b9b97703SBarry Smith                                         MatView_MPIDense,
1148357abbc8SBarry Smith                                         0,
114997304618SKris Buschelman                                         0,
1150d519adbfSMatthew Knepley                                 /* 64*/ 0,
115197304618SKris Buschelman                                         0,
115297304618SKris Buschelman                                         0,
115397304618SKris Buschelman                                         0,
115497304618SKris Buschelman                                         0,
1155d519adbfSMatthew Knepley                                 /* 69*/ 0,
115697304618SKris Buschelman                                         0,
115797304618SKris Buschelman                                         0,
115897304618SKris Buschelman                                         0,
115997304618SKris Buschelman                                         0,
1160d519adbfSMatthew Knepley                                 /* 74*/ 0,
116197304618SKris Buschelman                                         0,
116297304618SKris Buschelman                                         0,
116397304618SKris Buschelman                                         0,
116497304618SKris Buschelman                                         0,
1165d519adbfSMatthew Knepley                                 /* 79*/ 0,
116697304618SKris Buschelman                                         0,
116797304618SKris Buschelman                                         0,
116897304618SKris Buschelman                                         0,
11695bba2384SShri Abhyankar                                 /* 83*/ MatLoad_MPIDense,
1170865e5f61SKris Buschelman                                         0,
1171865e5f61SKris Buschelman                                         0,
1172865e5f61SKris Buschelman                                         0,
1173865e5f61SKris Buschelman                                         0,
1174865e5f61SKris Buschelman                                         0,
1175d519adbfSMatthew Knepley                                 /* 89*/
1176865e5f61SKris Buschelman                                         0,
1177865e5f61SKris Buschelman                                         0,
1178865e5f61SKris Buschelman                                         0,
11792fbe02b9SBarry Smith                                         0,
1180ba337c44SJed Brown                                         0,
1181d519adbfSMatthew Knepley                                 /* 94*/ 0,
1182865e5f61SKris Buschelman                                         0,
1183865e5f61SKris Buschelman                                         0,
1184ba337c44SJed Brown                                         0,
1185ba337c44SJed Brown                                         0,
1186ba337c44SJed Brown                                 /* 99*/ 0,
1187ba337c44SJed Brown                                         0,
1188ba337c44SJed Brown                                         0,
1189ba337c44SJed Brown                                         MatConjugate_MPIDense,
1190ba337c44SJed Brown                                         0,
1191ba337c44SJed Brown                                 /*104*/ 0,
1192ba337c44SJed Brown                                         MatRealPart_MPIDense,
1193ba337c44SJed Brown                                         MatImaginaryPart_MPIDense,
119486d161a7SShri Abhyankar                                         0,
119586d161a7SShri Abhyankar                                         0,
119686d161a7SShri Abhyankar                                 /*109*/ 0,
119786d161a7SShri Abhyankar                                         0,
119886d161a7SShri Abhyankar                                         0,
119986d161a7SShri Abhyankar                                         0,
120086d161a7SShri Abhyankar                                         0,
120186d161a7SShri Abhyankar                                 /*114*/ 0,
120286d161a7SShri Abhyankar                                         0,
120386d161a7SShri Abhyankar                                         0,
120486d161a7SShri Abhyankar                                         0,
120586d161a7SShri Abhyankar                                         0,
120686d161a7SShri Abhyankar                                 /*119*/ 0,
120786d161a7SShri Abhyankar                                         0,
120886d161a7SShri Abhyankar                                         0,
12090716a85fSBarry Smith                                         0,
12100716a85fSBarry Smith                                         0,
12110716a85fSBarry Smith                                 /*124*/ 0,
12120716a85fSBarry Smith                                         MatGetColumnNorms_MPIDense
1213ba337c44SJed Brown };
12148965ea79SLois Curfman McInnes 
1215273d9f13SBarry Smith EXTERN_C_BEGIN
12164a2ae208SSatish Balay #undef __FUNCT__
1217a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
12187087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1219a23d5eceSKris Buschelman {
1220a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1221dfbe8321SBarry Smith   PetscErrorCode ierr;
1222a23d5eceSKris Buschelman 
1223a23d5eceSKris Buschelman   PetscFunctionBegin;
1224a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1225a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1226a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1227a23d5eceSKris Buschelman 
1228a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
122934ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
123034ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
123134ef9618SShri Abhyankar   a->nvec = mat->cmap->n;
123234ef9618SShri Abhyankar 
1233f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1234d0f46423SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr);
12355c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
12365c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
123752e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1238a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1239a23d5eceSKris Buschelman }
1240a23d5eceSKris Buschelman EXTERN_C_END
1241a23d5eceSKris Buschelman 
1242a23d5eceSKris Buschelman EXTERN_C_BEGIN
1243a23d5eceSKris Buschelman #undef __FUNCT__
12444a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
12457087cfbeSBarry Smith PetscErrorCode  MatCreate_MPIDense(Mat mat)
1246273d9f13SBarry Smith {
1247273d9f13SBarry Smith   Mat_MPIDense   *a;
1248dfbe8321SBarry Smith   PetscErrorCode ierr;
1249273d9f13SBarry Smith 
1250273d9f13SBarry Smith   PetscFunctionBegin;
125138f2d2fdSLisandro Dalcin   ierr              = PetscNewLog(mat,Mat_MPIDense,&a);CHKERRQ(ierr);
1252b0a32e0cSBarry Smith   mat->data         = (void*)a;
1253273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1254273d9f13SBarry Smith 
1255273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
12567adad957SLisandro Dalcin   ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&a->rank);CHKERRQ(ierr);
12577adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)mat)->comm,&a->size);CHKERRQ(ierr);
1258273d9f13SBarry Smith 
1259273d9f13SBarry Smith   /* build cache for off array entries formed */
1260273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
12617adad957SLisandro Dalcin   ierr = MatStashCreate_Private(((PetscObject)mat)->comm,1,&mat->stash);CHKERRQ(ierr);
1262273d9f13SBarry Smith 
1263273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1264273d9f13SBarry Smith   a->lvec        = 0;
1265273d9f13SBarry Smith   a->Mvctx       = 0;
1266273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1267273d9f13SBarry Smith 
12688c778c55SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseGetArray_C","MatDenseGetArray_MPIDense",MatDenseGetArray_MPIDense);CHKERRQ(ierr);
12698c778c55SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseRestoreArray_C","MatDenseRestoreArray_MPIDense",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr);
12708c778c55SBarry Smith 
1271273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1272273d9f13SBarry Smith                                            "MatGetDiagonalBlock_MPIDense",
1273273d9f13SBarry Smith                                             MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1274a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1275a23d5eceSKris Buschelman                                            "MatMPIDenseSetPreallocation_MPIDense",
1276a23d5eceSKris Buschelman                                             MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
12774ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",
12784ae313f4SHong Zhang                                            "MatMatMult_MPIAIJ_MPIDense",
12794ae313f4SHong Zhang                                             MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
12804ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",
12814ae313f4SHong Zhang                                             "MatMatMultSymbolic_MPIAIJ_MPIDense",
12824ae313f4SHong Zhang                                              MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
12834ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",
12844ae313f4SHong Zhang                                            "MatMatMultNumeric_MPIAIJ_MPIDense",
12854ae313f4SHong Zhang                                             MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
128638aed534SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr);
128701b82886SBarry Smith 
1288273d9f13SBarry Smith   PetscFunctionReturn(0);
1289273d9f13SBarry Smith }
1290273d9f13SBarry Smith EXTERN_C_END
1291273d9f13SBarry Smith 
1292209238afSKris Buschelman /*MC
1293002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1294209238afSKris Buschelman 
1295209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1296209238afSKris Buschelman    and MATMPIDENSE otherwise.
1297209238afSKris Buschelman 
1298209238afSKris Buschelman    Options Database Keys:
1299209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1300209238afSKris Buschelman 
1301209238afSKris Buschelman   Level: beginner
1302209238afSKris Buschelman 
130301b82886SBarry Smith 
1304209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1305209238afSKris Buschelman M*/
1306209238afSKris Buschelman 
13074a2ae208SSatish Balay #undef __FUNCT__
13084a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1309273d9f13SBarry Smith /*@C
1310273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1311273d9f13SBarry Smith 
1312273d9f13SBarry Smith    Not collective
1313273d9f13SBarry Smith 
1314273d9f13SBarry Smith    Input Parameters:
1315273d9f13SBarry Smith .  A - the matrix
1316273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1317273d9f13SBarry Smith    to control all matrix memory allocation.
1318273d9f13SBarry Smith 
1319273d9f13SBarry Smith    Notes:
1320273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1321273d9f13SBarry Smith    storage by columns.
1322273d9f13SBarry Smith 
1323273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1324273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1325273d9f13SBarry Smith    set data=PETSC_NULL.
1326273d9f13SBarry Smith 
1327273d9f13SBarry Smith    Level: intermediate
1328273d9f13SBarry Smith 
1329273d9f13SBarry Smith .keywords: matrix,dense, parallel
1330273d9f13SBarry Smith 
1331273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1332273d9f13SBarry Smith @*/
13337087cfbeSBarry Smith PetscErrorCode  MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1334273d9f13SBarry Smith {
13354ac538c5SBarry Smith   PetscErrorCode ierr;
1336273d9f13SBarry Smith 
1337273d9f13SBarry Smith   PetscFunctionBegin;
13384ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar *),(mat,data));CHKERRQ(ierr);
1339273d9f13SBarry Smith   PetscFunctionReturn(0);
1340273d9f13SBarry Smith }
1341273d9f13SBarry Smith 
13424a2ae208SSatish Balay #undef __FUNCT__
134369b1f4b7SBarry Smith #define __FUNCT__ "MatCreateDense"
13448965ea79SLois Curfman McInnes /*@C
134569b1f4b7SBarry Smith    MatCreateDense - Creates a parallel matrix in dense format.
13468965ea79SLois Curfman McInnes 
1347db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1348db81eaa0SLois Curfman McInnes 
13498965ea79SLois Curfman McInnes    Input Parameters:
1350db81eaa0SLois Curfman McInnes +  comm - MPI communicator
13518965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1352db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
13538965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1354db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
13557f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1356dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
13578965ea79SLois Curfman McInnes 
13588965ea79SLois Curfman McInnes    Output Parameter:
1359477f1c0bSLois Curfman McInnes .  A - the matrix
13608965ea79SLois Curfman McInnes 
1361b259b22eSLois Curfman McInnes    Notes:
136239ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
136339ddd567SLois Curfman McInnes    storage by columns.
13648965ea79SLois Curfman McInnes 
136518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
136618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
13677f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
136818f449edSLois Curfman McInnes 
13698965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
13708965ea79SLois Curfman McInnes    (possibly both).
13718965ea79SLois Curfman McInnes 
1372027ccd11SLois Curfman McInnes    Level: intermediate
1373027ccd11SLois Curfman McInnes 
137439ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
13758965ea79SLois Curfman McInnes 
137639ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
13778965ea79SLois Curfman McInnes @*/
137869b1f4b7SBarry Smith PetscErrorCode  MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
13798965ea79SLois Curfman McInnes {
13806849ba73SBarry Smith   PetscErrorCode ierr;
138113f74950SBarry Smith   PetscMPIInt    size;
13828965ea79SLois Curfman McInnes 
13833a40ed3dSBarry Smith   PetscFunctionBegin;
1384f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1385f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1386273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1387273d9f13SBarry Smith   if (size > 1) {
1388273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1389273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1390273d9f13SBarry Smith   } else {
1391273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1392273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
13938c469469SLois Curfman McInnes   }
13943a40ed3dSBarry Smith   PetscFunctionReturn(0);
13958965ea79SLois Curfman McInnes }
13968965ea79SLois Curfman McInnes 
13974a2ae208SSatish Balay #undef __FUNCT__
13984a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
13996849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
14008965ea79SLois Curfman McInnes {
14018965ea79SLois Curfman McInnes   Mat            mat;
14023501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1403dfbe8321SBarry Smith   PetscErrorCode ierr;
14048965ea79SLois Curfman McInnes 
14053a40ed3dSBarry Smith   PetscFunctionBegin;
14068965ea79SLois Curfman McInnes   *newmat       = 0;
14077adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&mat);CHKERRQ(ierr);
1408d0f46423SBarry Smith   ierr = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
14097adad957SLisandro Dalcin   ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1410834f8fabSBarry Smith   a                 = (Mat_MPIDense*)mat->data;
1411e04c1aa4SHong Zhang   ierr              = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
14125aa7edbeSHong Zhang 
1413d5f3da31SBarry Smith   mat->factortype   = A->factortype;
1414c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1415273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
14168965ea79SLois Curfman McInnes 
14178965ea79SLois Curfman McInnes   a->size           = oldmat->size;
14188965ea79SLois Curfman McInnes   a->rank           = oldmat->rank;
1419e0fa3b82SLois Curfman McInnes   mat->insertmode   = NOT_SET_VALUES;
1420b9b97703SBarry Smith   a->nvec           = oldmat->nvec;
14213782ba37SSatish Balay   a->donotstash     = oldmat->donotstash;
1422e04c1aa4SHong Zhang 
14231e1e43feSBarry Smith   ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr);
14241e1e43feSBarry Smith   ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr);
14258965ea79SLois Curfman McInnes 
1426329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
14275609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
142852e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
142901b82886SBarry Smith 
14308965ea79SLois Curfman McInnes   *newmat = mat;
14313a40ed3dSBarry Smith   PetscFunctionReturn(0);
14328965ea79SLois Curfman McInnes }
14338965ea79SLois Curfman McInnes 
14344a2ae208SSatish Balay #undef __FUNCT__
14355bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
14365bba2384SShri Abhyankar PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat,PetscInt sizesset)
143786d161a7SShri Abhyankar {
143886d161a7SShri Abhyankar   PetscErrorCode ierr;
143986d161a7SShri Abhyankar   PetscMPIInt    rank,size;
144086d161a7SShri Abhyankar   PetscInt       *rowners,i,m,nz,j;
144186d161a7SShri Abhyankar   PetscScalar    *array,*vals,*vals_ptr;
144286d161a7SShri Abhyankar 
144386d161a7SShri Abhyankar   PetscFunctionBegin;
144486d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
144586d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
144686d161a7SShri Abhyankar 
144786d161a7SShri Abhyankar   /* determine ownership of all rows */
144886d161a7SShri Abhyankar   if (newmat->rmap->n < 0) m          = M/size + ((M % size) > rank);
144986d161a7SShri Abhyankar   else m = newmat->rmap->n;
145086d161a7SShri Abhyankar   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
145186d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
145286d161a7SShri Abhyankar   rowners[0] = 0;
145386d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
145486d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
145586d161a7SShri Abhyankar   }
145686d161a7SShri Abhyankar 
145786d161a7SShri Abhyankar   if (!sizesset) {
145886d161a7SShri Abhyankar     ierr = MatSetSizes(newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
145986d161a7SShri Abhyankar   }
146086d161a7SShri Abhyankar   ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
14618c778c55SBarry Smith   ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr);
146286d161a7SShri Abhyankar 
146386d161a7SShri Abhyankar   if (!rank) {
146486d161a7SShri Abhyankar     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
146586d161a7SShri Abhyankar 
146686d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
146786d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
146886d161a7SShri Abhyankar 
146986d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
147086d161a7SShri Abhyankar     vals_ptr = vals;
147186d161a7SShri Abhyankar     for (i=0; i<m; i++) {
147286d161a7SShri Abhyankar       for (j=0; j<N; j++) {
147386d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
147486d161a7SShri Abhyankar       }
147586d161a7SShri Abhyankar     }
147686d161a7SShri Abhyankar 
147786d161a7SShri Abhyankar     /* read in other processors and ship out */
147886d161a7SShri Abhyankar     for (i=1; i<size; i++) {
147986d161a7SShri Abhyankar       nz   = (rowners[i+1] - rowners[i])*N;
148086d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1481a25532f0SBarry Smith       ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
148286d161a7SShri Abhyankar     }
148386d161a7SShri Abhyankar   } else {
148486d161a7SShri Abhyankar     /* receive numeric values */
148586d161a7SShri Abhyankar     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
148686d161a7SShri Abhyankar 
148786d161a7SShri Abhyankar     /* receive message of values*/
1488a25532f0SBarry Smith     ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr);
148986d161a7SShri Abhyankar 
149086d161a7SShri Abhyankar     /* insert into matrix-by row (this is why cannot directly read into array */
149186d161a7SShri Abhyankar     vals_ptr = vals;
149286d161a7SShri Abhyankar     for (i=0; i<m; i++) {
149386d161a7SShri Abhyankar       for (j=0; j<N; j++) {
149486d161a7SShri Abhyankar         array[i + j*m] = *vals_ptr++;
149586d161a7SShri Abhyankar       }
149686d161a7SShri Abhyankar     }
149786d161a7SShri Abhyankar   }
14988c778c55SBarry Smith   ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr);
149986d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
150086d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
150186d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150286d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150386d161a7SShri Abhyankar   PetscFunctionReturn(0);
150486d161a7SShri Abhyankar }
150586d161a7SShri Abhyankar 
150686d161a7SShri Abhyankar #undef __FUNCT__
15075bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense"
1508112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer)
150986d161a7SShri Abhyankar {
151086d161a7SShri Abhyankar   PetscScalar    *vals,*svals;
151186d161a7SShri Abhyankar   MPI_Comm       comm = ((PetscObject)viewer)->comm;
151286d161a7SShri Abhyankar   MPI_Status     status;
151386d161a7SShri Abhyankar   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
151486d161a7SShri Abhyankar   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
151586d161a7SShri Abhyankar   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
151686d161a7SShri Abhyankar   PetscInt       i,nz,j,rstart,rend,sizesset=1,grows,gcols;
151786d161a7SShri Abhyankar   int            fd;
151886d161a7SShri Abhyankar   PetscErrorCode ierr;
151986d161a7SShri Abhyankar 
152086d161a7SShri Abhyankar   PetscFunctionBegin;
152186d161a7SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
152286d161a7SShri Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
152386d161a7SShri Abhyankar   if (!rank) {
152486d161a7SShri Abhyankar     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
152586d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
152686d161a7SShri Abhyankar     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
152786d161a7SShri Abhyankar   }
152886d161a7SShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) sizesset = 0;
152986d161a7SShri Abhyankar 
153086d161a7SShri Abhyankar   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
153186d161a7SShri Abhyankar   M = header[1]; N = header[2]; nz = header[3];
153286d161a7SShri Abhyankar 
153386d161a7SShri Abhyankar   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
153486d161a7SShri Abhyankar   if (sizesset && newmat->rmap->N < 0) newmat->rmap->N = M;
153586d161a7SShri Abhyankar   if (sizesset && newmat->cmap->N < 0) newmat->cmap->N = N;
153686d161a7SShri Abhyankar 
153786d161a7SShri Abhyankar   /* If global sizes are set, check if they are consistent with that given in the file */
153886d161a7SShri Abhyankar   if (sizesset) {
153986d161a7SShri Abhyankar     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
154086d161a7SShri Abhyankar   }
1541abd38a8fSBarry 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);
1542abd38a8fSBarry 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);
154386d161a7SShri Abhyankar 
154486d161a7SShri Abhyankar   /*
154586d161a7SShri Abhyankar        Handle case where matrix is stored on disk as a dense matrix
154686d161a7SShri Abhyankar   */
154786d161a7SShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
15485bba2384SShri Abhyankar     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat,sizesset);CHKERRQ(ierr);
154986d161a7SShri Abhyankar     PetscFunctionReturn(0);
155086d161a7SShri Abhyankar   }
155186d161a7SShri Abhyankar 
155286d161a7SShri Abhyankar   /* determine ownership of all rows */
155386d161a7SShri Abhyankar   if (newmat->rmap->n < 0) m          = PetscMPIIntCast(M/size + ((M % size) > rank));
155486d161a7SShri Abhyankar   else m = PetscMPIIntCast(newmat->rmap->n);
155586d161a7SShri Abhyankar   ierr       = PetscMalloc((size+2)*sizeof(PetscMPIInt),&rowners);CHKERRQ(ierr);
155686d161a7SShri Abhyankar   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
155786d161a7SShri Abhyankar   rowners[0] = 0;
155886d161a7SShri Abhyankar   for (i=2; i<=size; i++) {
155986d161a7SShri Abhyankar     rowners[i] += rowners[i-1];
156086d161a7SShri Abhyankar   }
156186d161a7SShri Abhyankar   rstart = rowners[rank];
156286d161a7SShri Abhyankar   rend   = rowners[rank+1];
156386d161a7SShri Abhyankar 
156486d161a7SShri Abhyankar   /* distribute row lengths to all processors */
156586d161a7SShri Abhyankar   ierr    = PetscMalloc2(rend-rstart,PetscInt,&ourlens,rend-rstart,PetscInt,&offlens);CHKERRQ(ierr);
156686d161a7SShri Abhyankar   if (!rank) {
156786d161a7SShri Abhyankar     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
156886d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
156986d161a7SShri Abhyankar     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
157086d161a7SShri Abhyankar     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
157186d161a7SShri Abhyankar     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
157286d161a7SShri Abhyankar     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
157386d161a7SShri Abhyankar   } else {
157486d161a7SShri Abhyankar     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
157586d161a7SShri Abhyankar   }
157686d161a7SShri Abhyankar 
157786d161a7SShri Abhyankar   if (!rank) {
157886d161a7SShri Abhyankar     /* calculate the number of nonzeros on each processor */
157986d161a7SShri Abhyankar     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
158086d161a7SShri Abhyankar     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
158186d161a7SShri Abhyankar     for (i=0; i<size; i++) {
158286d161a7SShri Abhyankar       for (j=rowners[i]; j< rowners[i+1]; j++) {
158386d161a7SShri Abhyankar         procsnz[i] += rowlengths[j];
158486d161a7SShri Abhyankar       }
158586d161a7SShri Abhyankar     }
158686d161a7SShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
158786d161a7SShri Abhyankar 
158886d161a7SShri Abhyankar     /* determine max buffer needed and allocate it */
158986d161a7SShri Abhyankar     maxnz = 0;
159086d161a7SShri Abhyankar     for (i=0; i<size; i++) {
159186d161a7SShri Abhyankar       maxnz = PetscMax(maxnz,procsnz[i]);
159286d161a7SShri Abhyankar     }
159386d161a7SShri Abhyankar     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
159486d161a7SShri Abhyankar 
159586d161a7SShri Abhyankar     /* read in my part of the matrix column indices  */
159686d161a7SShri Abhyankar     nz = procsnz[0];
159786d161a7SShri Abhyankar     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
159886d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
159986d161a7SShri Abhyankar 
160086d161a7SShri Abhyankar     /* read in every one elses and ship off */
160186d161a7SShri Abhyankar     for (i=1; i<size; i++) {
160286d161a7SShri Abhyankar       nz   = procsnz[i];
160386d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
160486d161a7SShri Abhyankar       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
160586d161a7SShri Abhyankar     }
160686d161a7SShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
160786d161a7SShri Abhyankar   } else {
160886d161a7SShri Abhyankar     /* determine buffer space needed for message */
160986d161a7SShri Abhyankar     nz = 0;
161086d161a7SShri Abhyankar     for (i=0; i<m; i++) {
161186d161a7SShri Abhyankar       nz += ourlens[i];
161286d161a7SShri Abhyankar     }
161386d161a7SShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
161486d161a7SShri Abhyankar 
161586d161a7SShri Abhyankar     /* receive message of column indices*/
161686d161a7SShri Abhyankar     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
161786d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
161886d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
161986d161a7SShri Abhyankar   }
162086d161a7SShri Abhyankar 
162186d161a7SShri Abhyankar   /* loop over local rows, determining number of off diagonal entries */
162286d161a7SShri Abhyankar   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
162386d161a7SShri Abhyankar   jj = 0;
162486d161a7SShri Abhyankar   for (i=0; i<m; i++) {
162586d161a7SShri Abhyankar     for (j=0; j<ourlens[i]; j++) {
162686d161a7SShri Abhyankar       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
162786d161a7SShri Abhyankar       jj++;
162886d161a7SShri Abhyankar     }
162986d161a7SShri Abhyankar   }
163086d161a7SShri Abhyankar 
163186d161a7SShri Abhyankar   /* create our matrix */
163286d161a7SShri Abhyankar   for (i=0; i<m; i++) {
163386d161a7SShri Abhyankar     ourlens[i] -= offlens[i];
163486d161a7SShri Abhyankar   }
163586d161a7SShri Abhyankar 
163686d161a7SShri Abhyankar   if (!sizesset) {
163786d161a7SShri Abhyankar     ierr = MatSetSizes(newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
163886d161a7SShri Abhyankar   }
163986d161a7SShri Abhyankar   ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
164086d161a7SShri Abhyankar   for (i=0; i<m; i++) {
164186d161a7SShri Abhyankar     ourlens[i] += offlens[i];
164286d161a7SShri Abhyankar   }
164386d161a7SShri Abhyankar 
164486d161a7SShri Abhyankar   if (!rank) {
164586d161a7SShri Abhyankar     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
164686d161a7SShri Abhyankar 
164786d161a7SShri Abhyankar     /* read in my part of the matrix numerical values  */
164886d161a7SShri Abhyankar     nz = procsnz[0];
164986d161a7SShri Abhyankar     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
165086d161a7SShri Abhyankar 
165186d161a7SShri Abhyankar     /* insert into matrix */
165286d161a7SShri Abhyankar     jj      = rstart;
165386d161a7SShri Abhyankar     smycols = mycols;
165486d161a7SShri Abhyankar     svals   = vals;
165586d161a7SShri Abhyankar     for (i=0; i<m; i++) {
165686d161a7SShri Abhyankar       ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
165786d161a7SShri Abhyankar       smycols += ourlens[i];
165886d161a7SShri Abhyankar       svals   += ourlens[i];
165986d161a7SShri Abhyankar       jj++;
166086d161a7SShri Abhyankar     }
166186d161a7SShri Abhyankar 
166286d161a7SShri Abhyankar     /* read in other processors and ship out */
166386d161a7SShri Abhyankar     for (i=1; i<size; i++) {
166486d161a7SShri Abhyankar       nz   = procsnz[i];
166586d161a7SShri Abhyankar       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
166686d161a7SShri Abhyankar       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr);
166786d161a7SShri Abhyankar     }
166886d161a7SShri Abhyankar     ierr = PetscFree(procsnz);CHKERRQ(ierr);
166986d161a7SShri Abhyankar   } else {
167086d161a7SShri Abhyankar     /* receive numeric values */
167186d161a7SShri Abhyankar     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
167286d161a7SShri Abhyankar 
167386d161a7SShri Abhyankar     /* receive message of values*/
167486d161a7SShri Abhyankar     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr);
167586d161a7SShri Abhyankar     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
167686d161a7SShri Abhyankar     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
167786d161a7SShri Abhyankar 
167886d161a7SShri Abhyankar     /* insert into matrix */
167986d161a7SShri Abhyankar     jj      = rstart;
168086d161a7SShri Abhyankar     smycols = mycols;
168186d161a7SShri Abhyankar     svals   = vals;
168286d161a7SShri Abhyankar     for (i=0; i<m; i++) {
168386d161a7SShri Abhyankar       ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
168486d161a7SShri Abhyankar       smycols += ourlens[i];
168586d161a7SShri Abhyankar       svals   += ourlens[i];
168686d161a7SShri Abhyankar       jj++;
168786d161a7SShri Abhyankar     }
168886d161a7SShri Abhyankar   }
168986d161a7SShri Abhyankar   ierr = PetscFree2(ourlens,offlens);CHKERRQ(ierr);
169086d161a7SShri Abhyankar   ierr = PetscFree(vals);CHKERRQ(ierr);
169186d161a7SShri Abhyankar   ierr = PetscFree(mycols);CHKERRQ(ierr);
169286d161a7SShri Abhyankar   ierr = PetscFree(rowners);CHKERRQ(ierr);
169386d161a7SShri Abhyankar 
169486d161a7SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169586d161a7SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169686d161a7SShri Abhyankar   PetscFunctionReturn(0);
169786d161a7SShri Abhyankar }
169886d161a7SShri Abhyankar 
169986d161a7SShri Abhyankar #undef __FUNCT__
17006e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense"
1701ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool  *flag)
17026e4ee0c6SHong Zhang {
17036e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
17046e4ee0c6SHong Zhang   Mat            a,b;
1705ace3abfcSBarry Smith   PetscBool      flg;
17066e4ee0c6SHong Zhang   PetscErrorCode ierr;
170790ace30eSBarry Smith 
17086e4ee0c6SHong Zhang   PetscFunctionBegin;
17096e4ee0c6SHong Zhang   a = matA->A;
17106e4ee0c6SHong Zhang   b = matB->A;
17116e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
1712c3aae356SJed Brown   ierr = MPI_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr);
17136e4ee0c6SHong Zhang   PetscFunctionReturn(0);
17146e4ee0c6SHong Zhang }
171590ace30eSBarry Smith 
1716