xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision f1af5d2ffeae1f5fc391a89939f4818e47770ae3)
1*f1af5d2fSBarry Smith /*$Id: mpidense.c,v 1.130 1999/10/24 14:02:12 bsmith Exp bsmith $*/
28965ea79SLois Curfman McInnes 
3ed3cc1f0SBarry Smith /*
4ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
5ed3cc1f0SBarry Smith */
6ed3cc1f0SBarry Smith 
770f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
8f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
98965ea79SLois Curfman McInnes 
100de54da6SSatish Balay EXTERN_C_BEGIN
110de54da6SSatish Balay #undef __FUNC__
120de54da6SSatish Balay #define __FUNC__ "MatGetDiagonalBlock_MPIDense"
130de54da6SSatish Balay int MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
140de54da6SSatish Balay {
150de54da6SSatish Balay   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
160de54da6SSatish Balay   int          m = mdn->m,rstart = mdn->rstart,rank,ierr;
170de54da6SSatish Balay   Scalar       *array;
180de54da6SSatish Balay   MPI_Comm     comm;
190de54da6SSatish Balay 
200de54da6SSatish Balay   PetscFunctionBegin;
210de54da6SSatish Balay   if (mdn->M != mdn->N) SETERRQ(PETSC_ERR_SUP,0,"Only square matrices supported.");
220de54da6SSatish Balay 
230de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
240de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
250de54da6SSatish Balay 
260de54da6SSatish Balay   ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr);
270de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
280de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
290de54da6SSatish Balay   ierr = MatCreateSeqDense(comm,m,m,array+m*rstart,B);CHKERRQ(ierr);
300de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
310de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
320de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
330de54da6SSatish Balay 
340de54da6SSatish Balay   *iscopy = PETSC_TRUE;
350de54da6SSatish Balay   PetscFunctionReturn(0);
360de54da6SSatish Balay }
370de54da6SSatish Balay EXTERN_C_END
380de54da6SSatish Balay 
397ef1d9bdSSatish Balay extern int MatSetUpMultiply_MPIDense(Mat);
407ef1d9bdSSatish Balay 
415615d1e5SSatish Balay #undef __FUNC__
425615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIDense"
438f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
448965ea79SLois Curfman McInnes {
4539b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense *) mat->data;
4639b7565bSBarry Smith   int          ierr, i, j, rstart = A->rstart, rend = A->rend, row;
4739b7565bSBarry Smith   int          roworiented = A->roworiented;
488965ea79SLois Curfman McInnes 
493a40ed3dSBarry Smith   PetscFunctionBegin;
508965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
515ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
52a8c6a408SBarry Smith     if (idxm[i] >= A->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
538965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
548965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
5539b7565bSBarry Smith       if (roworiented) {
5639b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
573a40ed3dSBarry Smith       } else {
588965ea79SLois Curfman McInnes         for ( j=0; j<n; j++ ) {
595ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
60a8c6a408SBarry Smith           if (idxn[j] >= A->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
6139b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
6239b7565bSBarry Smith         }
638965ea79SLois Curfman McInnes       }
643a40ed3dSBarry Smith     } else {
653782ba37SSatish Balay       if ( !A->donotstash) {
6639b7565bSBarry Smith         if (roworiented) {
678798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
68d36fbae8SSatish Balay         } else {
698798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
7039b7565bSBarry Smith         }
71b49de8d1SLois Curfman McInnes       }
72b49de8d1SLois Curfman McInnes     }
733782ba37SSatish Balay   }
743a40ed3dSBarry Smith   PetscFunctionReturn(0);
75b49de8d1SLois Curfman McInnes }
76b49de8d1SLois Curfman McInnes 
775615d1e5SSatish Balay #undef __FUNC__
785615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIDense"
798f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
80b49de8d1SLois Curfman McInnes {
81b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
82b49de8d1SLois Curfman McInnes   int          ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row;
83b49de8d1SLois Curfman McInnes 
843a40ed3dSBarry Smith   PetscFunctionBegin;
85b49de8d1SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
86a8c6a408SBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
87a8c6a408SBarry Smith     if (idxm[i] >= mdn->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
88b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
89b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
90b49de8d1SLois Curfman McInnes       for ( j=0; j<n; j++ ) {
91a8c6a408SBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
92a8c6a408SBarry Smith         if (idxn[j] >= mdn->N) {
93a8c6a408SBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
94a8c6a408SBarry Smith         }
95b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
96b49de8d1SLois Curfman McInnes       }
97a8c6a408SBarry Smith     } else {
98a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported");
998965ea79SLois Curfman McInnes     }
1008965ea79SLois Curfman McInnes   }
1013a40ed3dSBarry Smith   PetscFunctionReturn(0);
1028965ea79SLois Curfman McInnes }
1038965ea79SLois Curfman McInnes 
1045615d1e5SSatish Balay #undef __FUNC__
1055615d1e5SSatish Balay #define __FUNC__ "MatGetArray_MPIDense"
1068f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array)
107ff14e315SSatish Balay {
108ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
109ff14e315SSatish Balay   int          ierr;
110ff14e315SSatish Balay 
1113a40ed3dSBarry Smith   PetscFunctionBegin;
112ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1133a40ed3dSBarry Smith   PetscFunctionReturn(0);
114ff14e315SSatish Balay }
115ff14e315SSatish Balay 
1165615d1e5SSatish Balay #undef __FUNC__
117ca3fa75bSLois Curfman McInnes #define __FUNC__ "MatGetSubMatrix_MPIDense"
118ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
119ca3fa75bSLois Curfman McInnes {
120ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data, *newmatd;
121ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense *) mat->A->data;
12272d926a5SLois Curfman McInnes   int          i, j, ierr, *irow, *icol, rstart, rend, nrows, ncols, nlrows, nlcols, rank;
123ca3fa75bSLois Curfman McInnes   Scalar       *av, *bv, *v = lmat->v;
124ca3fa75bSLois Curfman McInnes   Mat          newmat;
125ca3fa75bSLois Curfman McInnes 
126ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
1277eba5e9cSLois Curfman McInnes   ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr);
128ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
129ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
130ca3fa75bSLois Curfman McInnes   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
131ca3fa75bSLois Curfman McInnes   ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr);
132ca3fa75bSLois Curfman McInnes 
133ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1347eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1357eba5e9cSLois Curfman McInnes      original matrix! */
136ca3fa75bSLois Curfman McInnes 
137ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1387eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
139ca3fa75bSLois Curfman McInnes 
140ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
141ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
1427eba5e9cSLois Curfman McInnes     /* SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); */
1437eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
144ca3fa75bSLois Curfman McInnes     newmat = *B;
145ca3fa75bSLois Curfman McInnes   } else {
146ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
147ca3fa75bSLois Curfman McInnes     ierr = MatCreateMPIDense(A->comm,nrows,ncols,PETSC_DECIDE,PETSC_DECIDE,PETSC_NULL,&newmat);CHKERRQ(ierr);
148ca3fa75bSLois Curfman McInnes   }
149ca3fa75bSLois Curfman McInnes 
150ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
151ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense *) newmat->data;
152ca3fa75bSLois Curfman McInnes   bv = ((Mat_SeqDense *)newmatd->A->data)->v;
153ca3fa75bSLois Curfman McInnes 
154ca3fa75bSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) {
155ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
156ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++ ) {
1577eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
158ca3fa75bSLois Curfman McInnes     }
159ca3fa75bSLois Curfman McInnes   }
160ca3fa75bSLois Curfman McInnes 
161ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
162ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
163ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164ca3fa75bSLois Curfman McInnes 
165ca3fa75bSLois Curfman McInnes   /* Free work space */
166ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
167ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
168ca3fa75bSLois Curfman McInnes   *B = newmat;
169ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
170ca3fa75bSLois Curfman McInnes }
171ca3fa75bSLois Curfman McInnes 
172ca3fa75bSLois Curfman McInnes #undef __FUNC__
1735615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_MPIDense"
1748f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array)
175ff14e315SSatish Balay {
1763a40ed3dSBarry Smith   PetscFunctionBegin;
1773a40ed3dSBarry Smith   PetscFunctionReturn(0);
178ff14e315SSatish Balay }
179ff14e315SSatish Balay 
1805615d1e5SSatish Balay #undef __FUNC__
1815615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIDense"
1828f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
1838965ea79SLois Curfman McInnes {
18439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
1858965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
186d36fbae8SSatish Balay   int          ierr,nstash,reallocs;
1878965ea79SLois Curfman McInnes   InsertMode   addv;
1888965ea79SLois Curfman McInnes 
1893a40ed3dSBarry Smith   PetscFunctionBegin;
1908965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
191ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
1927056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
193a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Cannot mix adds/inserts on different procs");
1948965ea79SLois Curfman McInnes   }
195e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
1968965ea79SLois Curfman McInnes 
1978798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
1988798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
199d36fbae8SSatish Balay   PLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",
200d36fbae8SSatish Balay            nstash,reallocs);
2013a40ed3dSBarry Smith   PetscFunctionReturn(0);
2028965ea79SLois Curfman McInnes }
2038965ea79SLois Curfman McInnes 
2045615d1e5SSatish Balay #undef __FUNC__
2055615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIDense"
2068f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2078965ea79SLois Curfman McInnes {
20839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
2097ef1d9bdSSatish Balay   int          i,n,ierr,*row,*col,flg,j,rstart,ncols;
2107ef1d9bdSSatish Balay   Scalar       *val;
211e0fa3b82SLois Curfman McInnes   InsertMode   addv=mat->insertmode;
2128965ea79SLois Curfman McInnes 
2133a40ed3dSBarry Smith   PetscFunctionBegin;
2148965ea79SLois Curfman McInnes   /*  wait on receives */
2157ef1d9bdSSatish Balay   while (1) {
2168798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2177ef1d9bdSSatish Balay     if (!flg) break;
2188965ea79SLois Curfman McInnes 
2197ef1d9bdSSatish Balay     for ( i=0; i<n; ) {
2207ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2217ef1d9bdSSatish Balay       for ( j=i,rstart=row[j]; j<n; j++ ) { if (row[j] != rstart) break; }
2227ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2237ef1d9bdSSatish Balay       else       ncols = n-i;
2247ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2257ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2267ef1d9bdSSatish Balay       i = j;
2278965ea79SLois Curfman McInnes     }
2287ef1d9bdSSatish Balay   }
2298798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2308965ea79SLois Curfman McInnes 
23139ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
23239ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2338965ea79SLois Curfman McInnes 
2346d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
23539ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2368965ea79SLois Curfman McInnes   }
2373a40ed3dSBarry Smith   PetscFunctionReturn(0);
2388965ea79SLois Curfman McInnes }
2398965ea79SLois Curfman McInnes 
2405615d1e5SSatish Balay #undef __FUNC__
2415615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIDense"
2428f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A)
2438965ea79SLois Curfman McInnes {
2443a40ed3dSBarry Smith   int          ierr;
24539ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense *) A->data;
2463a40ed3dSBarry Smith 
2473a40ed3dSBarry Smith   PetscFunctionBegin;
2483a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2493a40ed3dSBarry Smith   PetscFunctionReturn(0);
2508965ea79SLois Curfman McInnes }
2518965ea79SLois Curfman McInnes 
2525615d1e5SSatish Balay #undef __FUNC__
2535615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_MPIDense"
2548f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs)
2554e220ebcSLois Curfman McInnes {
2563a40ed3dSBarry Smith   PetscFunctionBegin;
2574e220ebcSLois Curfman McInnes   *bs = 1;
2583a40ed3dSBarry Smith   PetscFunctionReturn(0);
2594e220ebcSLois Curfman McInnes }
2604e220ebcSLois Curfman McInnes 
2618965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2628965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2638965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2643501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2658965ea79SLois Curfman McInnes    routine.
2668965ea79SLois Curfman McInnes */
2675615d1e5SSatish Balay #undef __FUNC__
2685615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIDense"
2698f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag)
2708965ea79SLois Curfman McInnes {
27139ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense *) A->data;
2728965ea79SLois Curfman McInnes   int            i,ierr,N, *rows,*owners = l->rowners,size = l->size;
2738965ea79SLois Curfman McInnes   int            *procs,*nprocs,j,found,idx,nsends,*work;
2748965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2758965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2768965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2778965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2788965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2798965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2808965ea79SLois Curfman McInnes   IS             istmp;
2818965ea79SLois Curfman McInnes 
2823a40ed3dSBarry Smith   PetscFunctionBegin;
28377c4ece6SBarry Smith   ierr = ISGetSize(is,&N);CHKERRQ(ierr);
2848965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
2858965ea79SLois Curfman McInnes 
2868965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
2870452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) );CHKPTRQ(nprocs);
288549d3d68SSatish Balay   ierr   = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
289549d3d68SSatish Balay   procs  = nprocs + size;
2900452661fSBarry Smith   owner  = (int *) PetscMalloc((N+1)*sizeof(int));CHKPTRQ(owner); /* see note*/
2918965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
2928965ea79SLois Curfman McInnes     idx = rows[i];
2938965ea79SLois Curfman McInnes     found = 0;
2948965ea79SLois Curfman McInnes     for ( j=0; j<size; j++ ) {
2958965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
2968965ea79SLois Curfman McInnes         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
2978965ea79SLois Curfman McInnes       }
2988965ea79SLois Curfman McInnes     }
299a8c6a408SBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range");
3008965ea79SLois Curfman McInnes   }
3018965ea79SLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
3028965ea79SLois Curfman McInnes 
3038965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
3046831982aSBarry Smith   work   = (int *) PetscMalloc( 2*size*sizeof(int) );CHKPTRQ(work);
3056831982aSBarry Smith   ierr   = MPI_Allreduce( nprocs, work,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr);
3068965ea79SLois Curfman McInnes   nmax   = work[rank];
3076831982aSBarry Smith   nrecvs = work[size+rank];
308606d414cSSatish Balay   ierr   = PetscFree(work);CHKERRQ(ierr);
3098965ea79SLois Curfman McInnes 
3108965ea79SLois Curfman McInnes   /* post receives:   */
3113a40ed3dSBarry Smith   rvalues    = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues);
3123a40ed3dSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits);
3138965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
314ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3158965ea79SLois Curfman McInnes   }
3168965ea79SLois Curfman McInnes 
3178965ea79SLois Curfman McInnes   /* do sends:
3188965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3198965ea79SLois Curfman McInnes          the ith processor
3208965ea79SLois Curfman McInnes   */
3210452661fSBarry Smith   svalues    = (int *) PetscMalloc( (N+1)*sizeof(int) );CHKPTRQ(svalues);
3227056b6fcSBarry Smith   send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
3230452661fSBarry Smith   starts     = (int *) PetscMalloc( (size+1)*sizeof(int) );CHKPTRQ(starts);
3248965ea79SLois Curfman McInnes   starts[0]  = 0;
3258965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3268965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
3278965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3288965ea79SLois Curfman McInnes   }
3298965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3308965ea79SLois Curfman McInnes 
3318965ea79SLois Curfman McInnes   starts[0] = 0;
3328965ea79SLois Curfman McInnes   for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3338965ea79SLois Curfman McInnes   count = 0;
3348965ea79SLois Curfman McInnes   for ( i=0; i<size; i++ ) {
3358965ea79SLois Curfman McInnes     if (procs[i]) {
336ca161407SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3378965ea79SLois Curfman McInnes     }
3388965ea79SLois Curfman McInnes   }
339606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3408965ea79SLois Curfman McInnes 
3418965ea79SLois Curfman McInnes   base = owners[rank];
3428965ea79SLois Curfman McInnes 
3438965ea79SLois Curfman McInnes   /*  wait on receives */
3440452661fSBarry Smith   lens   = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) );CHKPTRQ(lens);
3458965ea79SLois Curfman McInnes   source = lens + nrecvs;
3468965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3478965ea79SLois Curfman McInnes   while (count) {
348ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3498965ea79SLois Curfman McInnes     /* unpack receives into our local space */
350ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3518965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3528965ea79SLois Curfman McInnes     lens[imdex]  = n;
3538965ea79SLois Curfman McInnes     slen += n;
3548965ea79SLois Curfman McInnes     count--;
3558965ea79SLois Curfman McInnes   }
356606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3578965ea79SLois Curfman McInnes 
3588965ea79SLois Curfman McInnes   /* move the data into the send scatter */
3590452661fSBarry Smith   lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) );CHKPTRQ(lrows);
3608965ea79SLois Curfman McInnes   count = 0;
3618965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
3628965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3638965ea79SLois Curfman McInnes     for ( j=0; j<lens[i]; j++ ) {
3648965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3658965ea79SLois Curfman McInnes     }
3668965ea79SLois Curfman McInnes   }
367606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
368606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
369606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
370606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3718965ea79SLois Curfman McInnes 
3728965ea79SLois Curfman McInnes   /* actually zap the local rows */
373029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
3748965ea79SLois Curfman McInnes   PLogObjectParent(A,istmp);
375606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3768965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3788965ea79SLois Curfman McInnes 
3798965ea79SLois Curfman McInnes   /* wait on sends */
3808965ea79SLois Curfman McInnes   if (nsends) {
3817056b6fcSBarry Smith     send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
382ca161407SBarry Smith     ierr        = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
383606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
3848965ea79SLois Curfman McInnes   }
385606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
386606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
3878965ea79SLois Curfman McInnes 
3883a40ed3dSBarry Smith   PetscFunctionReturn(0);
3898965ea79SLois Curfman McInnes }
3908965ea79SLois Curfman McInnes 
3915615d1e5SSatish Balay #undef __FUNC__
3925615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIDense"
3938f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3948965ea79SLois Curfman McInnes {
39539ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
3968965ea79SLois Curfman McInnes   int          ierr;
397c456f294SBarry Smith 
3983a40ed3dSBarry Smith   PetscFunctionBegin;
39943a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40043a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40144cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4023a40ed3dSBarry Smith   PetscFunctionReturn(0);
4038965ea79SLois Curfman McInnes }
4048965ea79SLois Curfman McInnes 
4055615d1e5SSatish Balay #undef __FUNC__
4065615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIDense"
4078f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4088965ea79SLois Curfman McInnes {
40939ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4108965ea79SLois Curfman McInnes   int          ierr;
411c456f294SBarry Smith 
4123a40ed3dSBarry Smith   PetscFunctionBegin;
41343a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41443a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41544cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4163a40ed3dSBarry Smith   PetscFunctionReturn(0);
4178965ea79SLois Curfman McInnes }
4188965ea79SLois Curfman McInnes 
4195615d1e5SSatish Balay #undef __FUNC__
4205615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIDense"
4218f6be9afSLois Curfman McInnes int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy)
422096963f5SLois Curfman McInnes {
423096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
424096963f5SLois Curfman McInnes   int          ierr;
4253501a2bdSLois Curfman McInnes   Scalar       zero = 0.0;
426096963f5SLois Curfman McInnes 
4273a40ed3dSBarry Smith   PetscFunctionBegin;
4283501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
42944cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
430537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
431537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4323a40ed3dSBarry Smith   PetscFunctionReturn(0);
433096963f5SLois Curfman McInnes }
434096963f5SLois Curfman McInnes 
4355615d1e5SSatish Balay #undef __FUNC__
4365615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIDense"
4378f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
438096963f5SLois Curfman McInnes {
439096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
440096963f5SLois Curfman McInnes   int          ierr;
441096963f5SLois Curfman McInnes 
4423a40ed3dSBarry Smith   PetscFunctionBegin;
4433501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
44444cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
445537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
446537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4473a40ed3dSBarry Smith   PetscFunctionReturn(0);
448096963f5SLois Curfman McInnes }
449096963f5SLois Curfman McInnes 
4505615d1e5SSatish Balay #undef __FUNC__
4515615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIDense"
4528f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v)
4538965ea79SLois Curfman McInnes {
45439ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
455096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data;
45644cd7ae7SLois Curfman McInnes   int          ierr, len, i, n, m = a->m, radd;
45744cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
458ed3cc1f0SBarry Smith 
4593a40ed3dSBarry Smith   PetscFunctionBegin;
46044cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
461096963f5SLois Curfman McInnes   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
462096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
463a8c6a408SBarry Smith   if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
46444cd7ae7SLois Curfman McInnes   len = PetscMin(aloc->m,aloc->n);
4657ddc982cSLois Curfman McInnes   radd = a->rstart*m;
46644cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
467096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
468096963f5SLois Curfman McInnes   }
4699a8c540fSSatish Balay   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
4703a40ed3dSBarry Smith   PetscFunctionReturn(0);
4718965ea79SLois Curfman McInnes }
4728965ea79SLois Curfman McInnes 
4735615d1e5SSatish Balay #undef __FUNC__
4745615d1e5SSatish Balay #define __FUNC__ "MatDestroy_MPIDense"
475e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat)
4768965ea79SLois Curfman McInnes {
4773501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4788965ea79SLois Curfman McInnes   int          ierr;
479ed3cc1f0SBarry Smith 
4803a40ed3dSBarry Smith   PetscFunctionBegin;
48194d884c6SBarry Smith 
48294d884c6SBarry Smith   if (mat->mapping) {
48394d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr);
48494d884c6SBarry Smith   }
48594d884c6SBarry Smith   if (mat->bmapping) {
48694d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr);
48794d884c6SBarry Smith   }
488aa482453SBarry Smith #if defined(PETSC_USE_LOG)
489e1311b90SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N);
4908965ea79SLois Curfman McInnes #endif
4918798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
492606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
4933501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
4943501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4953501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
496622d7880SLois Curfman McInnes   if (mdn->factor) {
497606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
498606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
499606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
500606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
501622d7880SLois Curfman McInnes   }
502606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
50361b13de0SBarry Smith   if (mat->rmap) {
50461b13de0SBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
50561b13de0SBarry Smith   }
50661b13de0SBarry Smith   if (mat->cmap) {
50761b13de0SBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
50890f02eecSBarry Smith   }
5098965ea79SLois Curfman McInnes   PLogObjectDestroy(mat);
5100452661fSBarry Smith   PetscHeaderDestroy(mat);
5113a40ed3dSBarry Smith   PetscFunctionReturn(0);
5128965ea79SLois Curfman McInnes }
51339ddd567SLois Curfman McInnes 
5145615d1e5SSatish Balay #undef __FUNC__
5155615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary"
51639ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer)
5178965ea79SLois Curfman McInnes {
51839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
5198965ea79SLois Curfman McInnes   int          ierr;
5207056b6fcSBarry Smith 
5213a40ed3dSBarry Smith   PetscFunctionBegin;
52239ddd567SLois Curfman McInnes   if (mdn->size == 1) {
52339ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5248965ea79SLois Curfman McInnes   }
525a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported");
5263a40ed3dSBarry Smith   PetscFunctionReturn(0);
5278965ea79SLois Curfman McInnes }
5288965ea79SLois Curfman McInnes 
5295615d1e5SSatish Balay #undef __FUNC__
530*f1af5d2fSBarry Smith #define __FUNC__ "MatView_MPIDense_ASCIIorDraworSocket"
531*f1af5d2fSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,Viewer viewer)
5328965ea79SLois Curfman McInnes {
53339ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
53477ed5343SBarry Smith   int          ierr, format, size = mdn->size, rank = mdn->rank;
53519bcc07fSBarry Smith   ViewerType   vtype;
536*f1af5d2fSBarry Smith   PetscTruth   isascii,isdraw;
5378965ea79SLois Curfman McInnes 
5383a40ed3dSBarry Smith   PetscFunctionBegin;
539*f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
540*f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
541*f1af5d2fSBarry Smith   if (isascii) {
5423a40ed3dSBarry Smith     ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
543888f2ed8SSatish Balay     ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
544639f9d9dSBarry Smith     if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
5454e220ebcSLois Curfman McInnes       MatInfo info;
546888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
5476831982aSBarry Smith       ierr = ViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m,
5486831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
5496831982aSBarry Smith       ierr = ViewerFlush(viewer);CHKERRQ(ierr);
5503501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5513a40ed3dSBarry Smith       PetscFunctionReturn(0);
55296f6c058SBarry Smith     } else if (format == VIEWER_FORMAT_ASCII_INFO) {
5533a40ed3dSBarry Smith       PetscFunctionReturn(0);
5548965ea79SLois Curfman McInnes     }
555*f1af5d2fSBarry Smith   } else if (isdraw) {
556*f1af5d2fSBarry Smith     Draw       draw;
557*f1af5d2fSBarry Smith     PetscTruth isnull;
558*f1af5d2fSBarry Smith 
559*f1af5d2fSBarry Smith     ierr = ViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
560*f1af5d2fSBarry Smith     ierr = DrawIsNull(draw, &isnull);CHKERRQ(ierr);
561*f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
562*f1af5d2fSBarry Smith   }
56377ed5343SBarry Smith 
5648965ea79SLois Curfman McInnes   if (size == 1) {
56539ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5663a40ed3dSBarry Smith   } else {
5678965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5688965ea79SLois Curfman McInnes     Mat          A;
56939ddd567SLois Curfman McInnes     int          M = mdn->M, N = mdn->N,m,row,i, nz, *cols;
57039ddd567SLois Curfman McInnes     Scalar       *vals;
57139ddd567SLois Curfman McInnes     Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data;
5728965ea79SLois Curfman McInnes 
5738965ea79SLois Curfman McInnes     if (!rank) {
574*f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5753a40ed3dSBarry Smith     } else {
576*f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5778965ea79SLois Curfman McInnes     }
5788965ea79SLois Curfman McInnes     PLogObjectParent(mat,A);
5798965ea79SLois Curfman McInnes 
58039ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
58139ddd567SLois Curfman McInnes        but it's quick for now */
58239ddd567SLois Curfman McInnes     row = mdn->rstart; m = Amdn->m;
5838965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
58439ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
58539ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
58639ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
58739ddd567SLois Curfman McInnes       row++;
5888965ea79SLois Curfman McInnes     }
5898965ea79SLois Curfman McInnes 
5906d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5916d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5928965ea79SLois Curfman McInnes     if (!rank) {
5936831982aSBarry Smith       Viewer sviewer;
5946831982aSBarry Smith       ierr = ViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
5956831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5966831982aSBarry Smith       ierr = ViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
5978965ea79SLois Curfman McInnes     }
5986831982aSBarry Smith     ierr = ViewerFlush(viewer);CHKERRQ(ierr);
5998965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
6008965ea79SLois Curfman McInnes   }
6013a40ed3dSBarry Smith   PetscFunctionReturn(0);
6028965ea79SLois Curfman McInnes }
6038965ea79SLois Curfman McInnes 
6045615d1e5SSatish Balay #undef __FUNC__
6055615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense"
606e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer)
6078965ea79SLois Curfman McInnes {
60839ddd567SLois Curfman McInnes   int        ierr;
609*f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
6108965ea79SLois Curfman McInnes 
611433994e6SBarry Smith   PetscFunctionBegin;
6120f5bd95cSBarry Smith 
6136831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
6146831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
615*f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
616*f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
6170f5bd95cSBarry Smith 
618*f1af5d2fSBarry Smith   if (isascii || issocket || isdraw) {
619*f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6200f5bd95cSBarry Smith   } else if (isbinary) {
6213a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6225cd90555SBarry Smith   } else {
6230f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6248965ea79SLois Curfman McInnes   }
6253a40ed3dSBarry Smith   PetscFunctionReturn(0);
6268965ea79SLois Curfman McInnes }
6278965ea79SLois Curfman McInnes 
6285615d1e5SSatish Balay #undef __FUNC__
6295615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense"
6308f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6318965ea79SLois Curfman McInnes {
6323501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6333501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6344e220ebcSLois Curfman McInnes   int          ierr;
6354e220ebcSLois Curfman McInnes   double       isend[5], irecv[5];
6368965ea79SLois Curfman McInnes 
6373a40ed3dSBarry Smith   PetscFunctionBegin;
6384e220ebcSLois Curfman McInnes   info->rows_global    = (double)mat->M;
6394e220ebcSLois Curfman McInnes   info->columns_global = (double)mat->N;
6404e220ebcSLois Curfman McInnes   info->rows_local     = (double)mat->m;
6414e220ebcSLois Curfman McInnes   info->columns_local  = (double)mat->N;
6424e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6434e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6444e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6454e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6468965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6474e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6484e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6494e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6504e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6514e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6528965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
653f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
6544e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6554e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6564e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6574e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6584e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6598965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
660f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
6614e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6624e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6634e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6644e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6654e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6668965ea79SLois Curfman McInnes   }
6674e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6684e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6694e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6703a40ed3dSBarry Smith   PetscFunctionReturn(0);
6718965ea79SLois Curfman McInnes }
6728965ea79SLois Curfman McInnes 
6738c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*);
6748aaee692SLois Curfman McInnes    extern int MatLUFactorNumeric_MPIDense(Mat,Mat*);
6758aaee692SLois Curfman McInnes    extern int MatLUFactor_MPIDense(Mat,IS,IS,double);
6768aaee692SLois Curfman McInnes    extern int MatSolve_MPIDense(Mat,Vec,Vec);
6778c469469SLois Curfman McInnes    extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec);
6788aaee692SLois Curfman McInnes    extern int MatSolveTrans_MPIDense(Mat,Vec,Vec);
6798aaee692SLois Curfman McInnes    extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */
6808aaee692SLois Curfman McInnes 
6815615d1e5SSatish Balay #undef __FUNC__
6825615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense"
6838f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6848965ea79SLois Curfman McInnes {
68539ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
6868965ea79SLois Curfman McInnes 
6873a40ed3dSBarry Smith   PetscFunctionBegin;
6886d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6896d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
6904787f768SSatish Balay       op == MAT_NEW_NONZERO_LOCATION_ERR ||
6914787f768SSatish Balay       op == MAT_NEW_NONZERO_ALLOCATION_ERR ||
692219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
693219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED) {
694b1fbbac0SLois Curfman McInnes         MatSetOption(a->A,op);
695b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
696aeafbbfcSLois Curfman McInnes         a->roworiented = 1;
6978965ea79SLois Curfman McInnes         MatSetOption(a->A,op);
698b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
699219d9a1aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
7006d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
7016d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
702b51ba29fSSatish Balay              op == MAT_YES_NEW_DIAGONALS ||
703b51ba29fSSatish Balay              op == MAT_USE_HASH_TABLE) {
704981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
7053a40ed3dSBarry Smith   } else if (op == MAT_COLUMN_ORIENTED) {
7063a40ed3dSBarry Smith     a->roworiented = 0; MatSetOption(a->A,op);
7073782ba37SSatish Balay   } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) {
7083782ba37SSatish Balay     a->donotstash = 1;
7093a40ed3dSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS) {
7103a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
7113a40ed3dSBarry Smith   } else {
7123a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
7133a40ed3dSBarry Smith   }
7143a40ed3dSBarry Smith   PetscFunctionReturn(0);
7158965ea79SLois Curfman McInnes }
7168965ea79SLois Curfman McInnes 
7175615d1e5SSatish Balay #undef __FUNC__
7185615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense"
7198f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n)
7208965ea79SLois Curfman McInnes {
7213501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7223a40ed3dSBarry Smith 
7233a40ed3dSBarry Smith   PetscFunctionBegin;
724*f1af5d2fSBarry Smith   if (m) *m = mat->M;
725*f1af5d2fSBarry Smith   if (n) *n = mat->N;
7263a40ed3dSBarry Smith   PetscFunctionReturn(0);
7278965ea79SLois Curfman McInnes }
7288965ea79SLois Curfman McInnes 
7295615d1e5SSatish Balay #undef __FUNC__
7305615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense"
7318f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n)
7328965ea79SLois Curfman McInnes {
7333501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7343a40ed3dSBarry Smith 
7353a40ed3dSBarry Smith   PetscFunctionBegin;
7368965ea79SLois Curfman McInnes   *m = mat->m; *n = mat->N;
7373a40ed3dSBarry Smith   PetscFunctionReturn(0);
7388965ea79SLois Curfman McInnes }
7398965ea79SLois Curfman McInnes 
7405615d1e5SSatish Balay #undef __FUNC__
7415615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense"
7428f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
7438965ea79SLois Curfman McInnes {
7443501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7453a40ed3dSBarry Smith 
7463a40ed3dSBarry Smith   PetscFunctionBegin;
7478965ea79SLois Curfman McInnes   *m = mat->rstart; *n = mat->rend;
7483a40ed3dSBarry Smith   PetscFunctionReturn(0);
7498965ea79SLois Curfman McInnes }
7508965ea79SLois Curfman McInnes 
7515615d1e5SSatish Balay #undef __FUNC__
7525615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense"
7538f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
7548965ea79SLois Curfman McInnes {
7553501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7563a40ed3dSBarry Smith   int          lrow, rstart = mat->rstart, rend = mat->rend,ierr;
7578965ea79SLois Curfman McInnes 
7583a40ed3dSBarry Smith   PetscFunctionBegin;
759a8c6a408SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows")
7608965ea79SLois Curfman McInnes   lrow = row - rstart;
7613a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7623a40ed3dSBarry Smith   PetscFunctionReturn(0);
7638965ea79SLois Curfman McInnes }
7648965ea79SLois Curfman McInnes 
7655615d1e5SSatish Balay #undef __FUNC__
7665615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense"
7678f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
7688965ea79SLois Curfman McInnes {
769606d414cSSatish Balay   int ierr;
770606d414cSSatish Balay 
7713a40ed3dSBarry Smith   PetscFunctionBegin;
772606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
773606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7743a40ed3dSBarry Smith   PetscFunctionReturn(0);
7758965ea79SLois Curfman McInnes }
7768965ea79SLois Curfman McInnes 
7775615d1e5SSatish Balay #undef __FUNC__
7785b2fa520SLois Curfman McInnes #define __FUNC__ "MatDiagonalScale_MPIDense"
7795b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7805b2fa520SLois Curfman McInnes {
7815b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
7825b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data;
7835b2fa520SLois Curfman McInnes   Scalar       *l,*r,x,*v;
78472d926a5SLois Curfman McInnes   int          ierr,i,j,s2a,s3a,s2,s3,m=mat->m,n=mat->n;
7855b2fa520SLois Curfman McInnes 
7865b2fa520SLois Curfman McInnes   PetscFunctionBegin;
78772d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7885b2fa520SLois Curfman McInnes   if (ll) {
78972d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
79072d926a5SLois Curfman McInnes     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector non-conforming local size");
7915b2fa520SLois Curfman McInnes     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7925b2fa520SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
7935b2fa520SLois Curfman McInnes       x = l[i];
7945b2fa520SLois Curfman McInnes       v = mat->v + i;
7955b2fa520SLois Curfman McInnes       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
7965b2fa520SLois Curfman McInnes     }
7975b2fa520SLois Curfman McInnes     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
7985b2fa520SLois Curfman McInnes     PLogFlops(n*m);
7995b2fa520SLois Curfman McInnes   }
8005b2fa520SLois Curfman McInnes   if (rr) {
80172d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
80272d926a5SLois Curfman McInnes     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec non-conforming local size");
8035b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8045b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8055b2fa520SLois Curfman McInnes     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8065b2fa520SLois Curfman McInnes     for ( i=0; i<n; i++ ) {
8075b2fa520SLois Curfman McInnes       x = r[i];
8085b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8095b2fa520SLois Curfman McInnes       for ( j=0; j<m; j++ ) { (*v++) *= x;}
8105b2fa520SLois Curfman McInnes     }
81172d926a5SLois Curfman McInnes     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
8125b2fa520SLois Curfman McInnes     PLogFlops(n*m);
8135b2fa520SLois Curfman McInnes   }
8145b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8155b2fa520SLois Curfman McInnes }
8165b2fa520SLois Curfman McInnes 
8175b2fa520SLois Curfman McInnes #undef __FUNC__
8185615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense"
8198f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm)
820096963f5SLois Curfman McInnes {
8213501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
8223501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data;
8233501a2bdSLois Curfman McInnes   int          ierr, i, j;
8243501a2bdSLois Curfman McInnes   double       sum = 0.0;
8253501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
8263501a2bdSLois Curfman McInnes 
8273a40ed3dSBarry Smith   PetscFunctionBegin;
8283501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
8293501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm);CHKERRQ(ierr);
8303501a2bdSLois Curfman McInnes   } else {
8313501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
8323501a2bdSLois Curfman McInnes       for (i=0; i<mat->n*mat->m; i++ ) {
833aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
834e20fef11SSatish Balay         sum += PetscReal(PetscConj(*v)*(*v)); v++;
8353501a2bdSLois Curfman McInnes #else
8363501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
8373501a2bdSLois Curfman McInnes #endif
8383501a2bdSLois Curfman McInnes       }
839ca161407SBarry Smith       ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
8403501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
8413501a2bdSLois Curfman McInnes       PLogFlops(2*mat->n*mat->m);
8423a40ed3dSBarry Smith     } else if (type == NORM_1) {
8433501a2bdSLois Curfman McInnes       double *tmp, *tmp2;
8440452661fSBarry Smith       tmp  = (double *) PetscMalloc( 2*mdn->N*sizeof(double) );CHKPTRQ(tmp);
8453501a2bdSLois Curfman McInnes       tmp2 = tmp + mdn->N;
846549d3d68SSatish Balay       ierr = PetscMemzero(tmp,2*mdn->N*sizeof(double));CHKERRQ(ierr);
847096963f5SLois Curfman McInnes       *norm = 0.0;
8483501a2bdSLois Curfman McInnes       v = mat->v;
8493501a2bdSLois Curfman McInnes       for ( j=0; j<mat->n; j++ ) {
8503501a2bdSLois Curfman McInnes         for ( i=0; i<mat->m; i++ ) {
85167e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8523501a2bdSLois Curfman McInnes         }
8533501a2bdSLois Curfman McInnes       }
854ca161407SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
8553501a2bdSLois Curfman McInnes       for ( j=0; j<mdn->N; j++ ) {
8563501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
8573501a2bdSLois Curfman McInnes       }
858606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
8593501a2bdSLois Curfman McInnes       PLogFlops(mat->n*mat->m);
8603a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
8613501a2bdSLois Curfman McInnes       double ntemp;
8623501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
863ca161407SBarry Smith       ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
8643a40ed3dSBarry Smith     } else {
865a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
8663501a2bdSLois Curfman McInnes     }
8673501a2bdSLois Curfman McInnes   }
8683a40ed3dSBarry Smith   PetscFunctionReturn(0);
8693501a2bdSLois Curfman McInnes }
8703501a2bdSLois Curfman McInnes 
8715615d1e5SSatish Balay #undef __FUNC__
8725615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense"
8738f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8743501a2bdSLois Curfman McInnes {
8753501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
8763501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data;
8773501a2bdSLois Curfman McInnes   Mat          B;
8783501a2bdSLois Curfman McInnes   int          M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart;
8793501a2bdSLois Curfman McInnes   int          j, i, ierr;
8803501a2bdSLois Curfman McInnes   Scalar       *v;
8813501a2bdSLois Curfman McInnes 
8823a40ed3dSBarry Smith   PetscFunctionBegin;
8837056b6fcSBarry Smith   if (matout == PETSC_NULL && M != N) {
884a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place");
8857056b6fcSBarry Smith   }
8867056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8873501a2bdSLois Curfman McInnes 
8883501a2bdSLois Curfman McInnes   m = Aloc->m; n = Aloc->n; v = Aloc->v;
8890452661fSBarry Smith   rwork = (int *) PetscMalloc(n*sizeof(int));CHKPTRQ(rwork);
8903501a2bdSLois Curfman McInnes   for ( j=0; j<n; j++ ) {
8913501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8923501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8933501a2bdSLois Curfman McInnes     v   += m;
8943501a2bdSLois Curfman McInnes   }
895606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8966d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8976d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8983638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
8993501a2bdSLois Curfman McInnes     *matout = B;
9003501a2bdSLois Curfman McInnes   } else {
901f830108cSBarry Smith     PetscOps *Abops;
90209dc0095SBarry Smith     MatOps   Aops;
903f830108cSBarry Smith 
9043501a2bdSLois Curfman McInnes     /* This isn't really an in-place transpose, but free data struct from a */
905606d414cSSatish Balay     ierr = PetscFree(a->rowners);CHKERRQ(ierr);
9063501a2bdSLois Curfman McInnes     ierr = MatDestroy(a->A);CHKERRQ(ierr);
9073501a2bdSLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
9083501a2bdSLois Curfman McInnes     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
909606d414cSSatish Balay     ierr = PetscFree(a);CHKERRQ(ierr);
910f830108cSBarry Smith 
911f830108cSBarry Smith     /*
912f830108cSBarry Smith          This is horrible, horrible code. We need to keep the
913f830108cSBarry Smith       A pointers for the bops and ops but copy everything
914f830108cSBarry Smith       else from C.
915f830108cSBarry Smith     */
916f830108cSBarry Smith     Abops   = A->bops;
917f830108cSBarry Smith     Aops    = A->ops;
918549d3d68SSatish Balay     ierr    = PetscMemcpy(A,B,sizeof(struct _p_Mat));CHKERRQ(ierr);
919f830108cSBarry Smith     A->bops = Abops;
920f830108cSBarry Smith     A->ops  = Aops;
921f830108cSBarry Smith 
9220452661fSBarry Smith     PetscHeaderDestroy(B);
9233501a2bdSLois Curfman McInnes   }
9243a40ed3dSBarry Smith   PetscFunctionReturn(0);
925096963f5SLois Curfman McInnes }
926096963f5SLois Curfman McInnes 
927eadb2fb4SBarry Smith #include "pinclude/blaslapack.h"
9285615d1e5SSatish Balay #undef __FUNC__
9295615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense"
9308f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA)
93144cd7ae7SLois Curfman McInnes {
93244cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense *) inA->data;
93344cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->A->data;
93444cd7ae7SLois Curfman McInnes   int          one = 1, nz;
93544cd7ae7SLois Curfman McInnes 
9363a40ed3dSBarry Smith   PetscFunctionBegin;
93744cd7ae7SLois Curfman McInnes   nz = a->m*a->n;
93844cd7ae7SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
93944cd7ae7SLois Curfman McInnes   PLogFlops(nz);
9403a40ed3dSBarry Smith   PetscFunctionReturn(0);
94144cd7ae7SLois Curfman McInnes }
94244cd7ae7SLois Curfman McInnes 
9435609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
9447b2a1423SBarry Smith extern int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **);
9458965ea79SLois Curfman McInnes 
9468965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
94709dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
94809dc0095SBarry Smith        MatGetRow_MPIDense,
94909dc0095SBarry Smith        MatRestoreRow_MPIDense,
95009dc0095SBarry Smith        MatMult_MPIDense,
95109dc0095SBarry Smith        MatMultAdd_MPIDense,
95209dc0095SBarry Smith        MatMultTrans_MPIDense,
95309dc0095SBarry Smith        MatMultTransAdd_MPIDense,
9548965ea79SLois Curfman McInnes        0,
95509dc0095SBarry Smith        0,
95609dc0095SBarry Smith        0,
95709dc0095SBarry Smith        0,
95809dc0095SBarry Smith        0,
95909dc0095SBarry Smith        0,
96009dc0095SBarry Smith        0,
96109dc0095SBarry Smith        MatTranspose_MPIDense,
96209dc0095SBarry Smith        MatGetInfo_MPIDense,0,
96309dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9645b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
96509dc0095SBarry Smith        MatNorm_MPIDense,
96609dc0095SBarry Smith        MatAssemblyBegin_MPIDense,
96709dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
96809dc0095SBarry Smith        0,
96909dc0095SBarry Smith        MatSetOption_MPIDense,
97009dc0095SBarry Smith        MatZeroEntries_MPIDense,
97109dc0095SBarry Smith        MatZeroRows_MPIDense,
97209dc0095SBarry Smith        0,
97309dc0095SBarry Smith        0,
97409dc0095SBarry Smith        0,
97509dc0095SBarry Smith        0,
97609dc0095SBarry Smith        MatGetSize_MPIDense,
97709dc0095SBarry Smith        MatGetLocalSize_MPIDense,
97839ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
97909dc0095SBarry Smith        0,
98009dc0095SBarry Smith        0,
98109dc0095SBarry Smith        MatGetArray_MPIDense,
98209dc0095SBarry Smith        MatRestoreArray_MPIDense,
9835609ef8eSBarry Smith        MatDuplicate_MPIDense,
98409dc0095SBarry Smith        0,
98509dc0095SBarry Smith        0,
98609dc0095SBarry Smith        0,
98709dc0095SBarry Smith        0,
98809dc0095SBarry Smith        0,
9892ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
99009dc0095SBarry Smith        0,
99109dc0095SBarry Smith        MatGetValues_MPIDense,
99209dc0095SBarry Smith        0,
99309dc0095SBarry Smith        0,
99409dc0095SBarry Smith        MatScale_MPIDense,
99509dc0095SBarry Smith        0,
99609dc0095SBarry Smith        0,
99709dc0095SBarry Smith        0,
99809dc0095SBarry Smith        MatGetBlockSize_MPIDense,
99909dc0095SBarry Smith        0,
100009dc0095SBarry Smith        0,
100109dc0095SBarry Smith        0,
100209dc0095SBarry Smith        0,
100309dc0095SBarry Smith        0,
100409dc0095SBarry Smith        0,
100509dc0095SBarry Smith        0,
100609dc0095SBarry Smith        0,
100709dc0095SBarry Smith        0,
1008ca3fa75bSLois Curfman McInnes        MatGetSubMatrix_MPIDense,
100909dc0095SBarry Smith        0,
101009dc0095SBarry Smith        0,
101109dc0095SBarry Smith        MatGetMaps_Petsc};
10128965ea79SLois Curfman McInnes 
10135615d1e5SSatish Balay #undef __FUNC__
10145615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense"
10158965ea79SLois Curfman McInnes /*@C
101639ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
10178965ea79SLois Curfman McInnes 
1018db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1019db81eaa0SLois Curfman McInnes 
10208965ea79SLois Curfman McInnes    Input Parameters:
1021db81eaa0SLois Curfman McInnes +  comm - MPI communicator
10228965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1023db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
10248965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1025db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1026db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1027dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
10288965ea79SLois Curfman McInnes 
10298965ea79SLois Curfman McInnes    Output Parameter:
1030477f1c0bSLois Curfman McInnes .  A - the matrix
10318965ea79SLois Curfman McInnes 
1032b259b22eSLois Curfman McInnes    Notes:
103339ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
103439ddd567SLois Curfman McInnes    storage by columns.
10358965ea79SLois Curfman McInnes 
103618f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
103718f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1038b4fd4287SBarry Smith    set data=PETSC_NULL.
103918f449edSLois Curfman McInnes 
10408965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
10418965ea79SLois Curfman McInnes    (possibly both).
10428965ea79SLois Curfman McInnes 
1043027ccd11SLois Curfman McInnes    Level: intermediate
1044027ccd11SLois Curfman McInnes 
104539ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel
10468965ea79SLois Curfman McInnes 
104739ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
10488965ea79SLois Curfman McInnes @*/
1049477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
10508965ea79SLois Curfman McInnes {
10518965ea79SLois Curfman McInnes   Mat          mat;
105239ddd567SLois Curfman McInnes   Mat_MPIDense *a;
1053*f1af5d2fSBarry Smith   int          ierr, i;
1054*f1af5d2fSBarry Smith   PetscTruth   flg;
10558965ea79SLois Curfman McInnes 
10563a40ed3dSBarry Smith   PetscFunctionBegin;
1057ed2daf61SLois Curfman McInnes   /* Note:  For now, when data is specified above, this assumes the user correctly
1058ed2daf61SLois Curfman McInnes    allocates the local dense storage space.  We should add error checking. */
105918f449edSLois Curfman McInnes 
1060477f1c0bSLois Curfman McInnes   *A = 0;
10613f1db9ecSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",comm,MatDestroy,MatView);
10628965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
10630452661fSBarry Smith   mat->data         = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a);
1064549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1065e1311b90SBarry Smith   mat->ops->destroy = MatDestroy_MPIDense;
1066e1311b90SBarry Smith   mat->ops->view    = MatView_MPIDense;
10678965ea79SLois Curfman McInnes   mat->factor       = 0;
106890f02eecSBarry Smith   mat->mapping      = 0;
10698965ea79SLois Curfman McInnes 
1070622d7880SLois Curfman McInnes   a->factor       = 0;
1071e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1072d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&a->rank);CHKERRQ(ierr);
1073d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&a->size);CHKERRQ(ierr);
10748965ea79SLois Curfman McInnes 
107596f6c058SBarry Smith   ierr = PetscSplitOwnership(comm,&m,&M);CHKERRQ(ierr);
107639ddd567SLois Curfman McInnes 
1077*f1af5d2fSBarry Smith   ierr = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr);
1078be0abb6dSBarry Smith   a->nvec = n;
1079c7fcc2eaSBarry Smith 
108039ddd567SLois Curfman McInnes   /* each row stores all columns */
1081aca0ad90SLois Curfman McInnes   a->N = mat->N = N;
1082aca0ad90SLois Curfman McInnes   a->M = mat->M = M;
1083aca0ad90SLois Curfman McInnes   a->m = mat->m = m;
1084be0abb6dSBarry Smith   a->n = mat->n = N;   /* NOTE: n == N */
10858965ea79SLois Curfman McInnes 
1086c7fcc2eaSBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1087c7fcc2eaSBarry Smith      we should remove the duplicate information that is not contained in the maps */
1088488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,M,&mat->rmap);CHKERRQ(ierr);
1089be0abb6dSBarry Smith   ierr = MapCreateMPI(comm,n,N,&mat->cmap);CHKERRQ(ierr);
1090c7fcc2eaSBarry Smith 
10918965ea79SLois Curfman McInnes   /* build local table of row and column ownerships */
1092d7e8b826SBarry Smith   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners);
1093d7e8b826SBarry Smith   a->cowners = a->rowners + a->size + 1;
1094f09e8eb9SSatish Balay   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1095ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
10968965ea79SLois Curfman McInnes   a->rowners[0] = 0;
10978965ea79SLois Curfman McInnes   for ( i=2; i<=a->size; i++ ) {
10988965ea79SLois Curfman McInnes     a->rowners[i] += a->rowners[i-1];
10998965ea79SLois Curfman McInnes   }
11008965ea79SLois Curfman McInnes   a->rstart = a->rowners[a->rank];
11018965ea79SLois Curfman McInnes   a->rend   = a->rowners[a->rank+1];
1102ca161407SBarry Smith   ierr      = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
1103d7e8b826SBarry Smith   a->cowners[0] = 0;
1104d7e8b826SBarry Smith   for ( i=2; i<=a->size; i++ ) {
1105d7e8b826SBarry Smith     a->cowners[i] += a->cowners[i-1];
1106d7e8b826SBarry Smith   }
11078965ea79SLois Curfman McInnes 
1108029af93fSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A);CHKERRQ(ierr);
11098965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
11108965ea79SLois Curfman McInnes 
11118965ea79SLois Curfman McInnes   /* build cache for off array entries formed */
11123782ba37SSatish Balay   a->donotstash = 0;
11138798bf22SSatish Balay   ierr = MatStashCreate_Private(comm,1,&mat->stash);CHKERRQ(ierr);
11148965ea79SLois Curfman McInnes 
11158965ea79SLois Curfman McInnes   /* stuff used for matrix vector multiply */
11168965ea79SLois Curfman McInnes   a->lvec        = 0;
11178965ea79SLois Curfman McInnes   a->Mvctx       = 0;
111839b7565bSBarry Smith   a->roworiented = 1;
11198965ea79SLois Curfman McInnes 
1120*f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
11210de54da6SSatish Balay                                      "MatGetDiagonalBlock_MPIDense",
11220de54da6SSatish Balay                                      (void*)MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
11230de54da6SSatish Balay 
1124477f1c0bSLois Curfman McInnes   *A = mat;
112525cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
112625cdf11fSBarry Smith   if (flg) {
11278c469469SLois Curfman McInnes     ierr = MatPrintHelp(mat);CHKERRQ(ierr);
11288c469469SLois Curfman McInnes   }
11293a40ed3dSBarry Smith   PetscFunctionReturn(0);
11308965ea79SLois Curfman McInnes }
11318965ea79SLois Curfman McInnes 
11325615d1e5SSatish Balay #undef __FUNC__
11335609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_MPIDense"
11345609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11358965ea79SLois Curfman McInnes {
11368965ea79SLois Curfman McInnes   Mat          mat;
11373501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data;
113839ddd567SLois Curfman McInnes   int          ierr;
11392ba99913SLois Curfman McInnes   FactorCtx    *factor;
11408965ea79SLois Curfman McInnes 
11413a40ed3dSBarry Smith   PetscFunctionBegin;
11428965ea79SLois Curfman McInnes   *newmat       = 0;
11433f1db9ecSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",A->comm,MatDestroy,MatView);
11448965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
11450452661fSBarry Smith   mat->data         = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a);
1146549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1147e1311b90SBarry Smith   mat->ops->destroy = MatDestroy_MPIDense;
1148e1311b90SBarry Smith   mat->ops->view    = MatView_MPIDense;
11493501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1150c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
11518965ea79SLois Curfman McInnes 
115244cd7ae7SLois Curfman McInnes   a->m = mat->m = oldmat->m;
115344cd7ae7SLois Curfman McInnes   a->n = mat->n = oldmat->n;
115444cd7ae7SLois Curfman McInnes   a->M = mat->M = oldmat->M;
115544cd7ae7SLois Curfman McInnes   a->N = mat->N = oldmat->N;
11562ba99913SLois Curfman McInnes   if (oldmat->factor) {
11572ba99913SLois Curfman McInnes     a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx));CHKPTRQ(factor);
11582ba99913SLois Curfman McInnes     /* copy factor contents ... add this code! */
11592ba99913SLois Curfman McInnes   } else a->factor = 0;
11608965ea79SLois Curfman McInnes 
11618965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
11628965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
11638965ea79SLois Curfman McInnes   a->size         = oldmat->size;
11648965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1165e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
11663782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
11670452661fSBarry Smith   a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners);
1168f09e8eb9SSatish Balay   PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1169549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
11708798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
11718965ea79SLois Curfman McInnes 
11728965ea79SLois Curfman McInnes   ierr =  VecDuplicate(oldmat->lvec,&a->lvec);CHKERRQ(ierr);
11738965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->lvec);
117455659b69SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx);CHKERRQ(ierr);
11758965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->Mvctx);
11765609ef8eSBarry Smith   ierr =  MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
11778965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
11788965ea79SLois Curfman McInnes   *newmat = mat;
11793a40ed3dSBarry Smith   PetscFunctionReturn(0);
11808965ea79SLois Curfman McInnes }
11818965ea79SLois Curfman McInnes 
118277c4ece6SBarry Smith #include "sys.h"
11838965ea79SLois Curfman McInnes 
11845615d1e5SSatish Balay #undef __FUNC__
11855615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile"
118690ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat)
118790ace30eSBarry Smith {
118840011551SBarry Smith   int        *rowners, i,size,rank,m,ierr,nz,j;
118990ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
119090ace30eSBarry Smith   MPI_Status status;
119190ace30eSBarry Smith 
11923a40ed3dSBarry Smith   PetscFunctionBegin;
1193d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1194d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
119590ace30eSBarry Smith 
119690ace30eSBarry Smith   /* determine ownership of all rows */
119790ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
119890ace30eSBarry Smith   rowners    = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners);
1199ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
120090ace30eSBarry Smith   rowners[0] = 0;
120190ace30eSBarry Smith   for ( i=2; i<=size; i++ ) {
120290ace30eSBarry Smith     rowners[i] += rowners[i-1];
120390ace30eSBarry Smith   }
120490ace30eSBarry Smith 
120590ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
120690ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
120790ace30eSBarry Smith 
120890ace30eSBarry Smith   if (!rank) {
120990ace30eSBarry Smith     vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals);
121090ace30eSBarry Smith 
121190ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
12120752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
121390ace30eSBarry Smith 
121490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
121590ace30eSBarry Smith     vals_ptr = vals;
121690ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
121790ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
121890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
121990ace30eSBarry Smith       }
122090ace30eSBarry Smith     }
122190ace30eSBarry Smith 
122290ace30eSBarry Smith     /* read in other processors and ship out */
122390ace30eSBarry Smith     for ( i=1; i<size; i++ ) {
122490ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12250752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1226ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
122790ace30eSBarry Smith     }
12283a40ed3dSBarry Smith   } else {
122990ace30eSBarry Smith     /* receive numeric values */
123090ace30eSBarry Smith     vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals);
123190ace30eSBarry Smith 
123290ace30eSBarry Smith     /* receive message of values*/
1233ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
123490ace30eSBarry Smith 
123590ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
123690ace30eSBarry Smith     vals_ptr = vals;
123790ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
123890ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
123990ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
124090ace30eSBarry Smith       }
124190ace30eSBarry Smith     }
124290ace30eSBarry Smith   }
1243606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1244606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
12456d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12466d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12473a40ed3dSBarry Smith   PetscFunctionReturn(0);
124890ace30eSBarry Smith }
124990ace30eSBarry Smith 
125090ace30eSBarry Smith 
12515615d1e5SSatish Balay #undef __FUNC__
12525615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense"
125319bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat)
12548965ea79SLois Curfman McInnes {
12558965ea79SLois Curfman McInnes   Mat          A;
12568965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
125719bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
12588965ea79SLois Curfman McInnes   MPI_Status   status;
12598965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
12608965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
126119bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
12623a40ed3dSBarry Smith   int          i, nz, ierr, j,rstart, rend, fd;
12638965ea79SLois Curfman McInnes 
12643a40ed3dSBarry Smith   PetscFunctionBegin;
1265d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1266d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12678965ea79SLois Curfman McInnes   if (!rank) {
126890ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
12690752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1270a8c6a408SBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object");
12718965ea79SLois Curfman McInnes   }
12728965ea79SLois Curfman McInnes 
1273ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
127490ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
127590ace30eSBarry Smith 
127690ace30eSBarry Smith   /*
127790ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
127890ace30eSBarry Smith   */
127990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
12803a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
12813a40ed3dSBarry Smith     PetscFunctionReturn(0);
128290ace30eSBarry Smith   }
128390ace30eSBarry Smith 
12848965ea79SLois Curfman McInnes   /* determine ownership of all rows */
12858965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
12860452661fSBarry Smith   rowners    = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners);
1287ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
12888965ea79SLois Curfman McInnes   rowners[0] = 0;
12898965ea79SLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
12908965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
12918965ea79SLois Curfman McInnes   }
12928965ea79SLois Curfman McInnes   rstart = rowners[rank];
12938965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
12948965ea79SLois Curfman McInnes 
12958965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
12960452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) );CHKPTRQ(ourlens);
12978965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
12988965ea79SLois Curfman McInnes   if (!rank) {
12990452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) );CHKPTRQ(rowlengths);
13000752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
13010452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(sndcounts);
13028965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
1303ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1304606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1305ca161407SBarry Smith   } else {
1306ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr);
13078965ea79SLois Curfman McInnes   }
13088965ea79SLois Curfman McInnes 
13098965ea79SLois Curfman McInnes   if (!rank) {
13108965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
13110452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(procsnz);
1312549d3d68SSatish Balay     ierr    = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
13138965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
13148965ea79SLois Curfman McInnes       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
13158965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
13168965ea79SLois Curfman McInnes       }
13178965ea79SLois Curfman McInnes     }
1318606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
13198965ea79SLois Curfman McInnes 
13208965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
13218965ea79SLois Curfman McInnes     maxnz = 0;
13228965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
13230452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
13248965ea79SLois Curfman McInnes     }
13250452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) );CHKPTRQ(cols);
13268965ea79SLois Curfman McInnes 
13278965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13288965ea79SLois Curfman McInnes     nz = procsnz[0];
13290452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols);
13300752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13318965ea79SLois Curfman McInnes 
13328965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13338965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
13348965ea79SLois Curfman McInnes       nz   = procsnz[i];
13350752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1336ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13378965ea79SLois Curfman McInnes     }
1338606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13393a40ed3dSBarry Smith   } else {
13408965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13418965ea79SLois Curfman McInnes     nz = 0;
13428965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
13438965ea79SLois Curfman McInnes       nz += ourlens[i];
13448965ea79SLois Curfman McInnes     }
13450452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols);
13468965ea79SLois Curfman McInnes 
13478965ea79SLois Curfman McInnes     /* receive message of column indices*/
1348ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1349ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
1350a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
13518965ea79SLois Curfman McInnes   }
13528965ea79SLois Curfman McInnes 
13538965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1354549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
13558965ea79SLois Curfman McInnes   jj = 0;
13568965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
13578965ea79SLois Curfman McInnes     for ( j=0; j<ourlens[i]; j++ ) {
13588965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
13598965ea79SLois Curfman McInnes       jj++;
13608965ea79SLois Curfman McInnes     }
13618965ea79SLois Curfman McInnes   }
13628965ea79SLois Curfman McInnes 
13638965ea79SLois Curfman McInnes   /* create our matrix */
13648965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
13658965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
13668965ea79SLois Curfman McInnes   }
1367b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
13688965ea79SLois Curfman McInnes   A = *newmat;
13698965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
13708965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
13718965ea79SLois Curfman McInnes   }
13728965ea79SLois Curfman McInnes 
13738965ea79SLois Curfman McInnes   if (!rank) {
13740452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) );CHKPTRQ(vals);
13758965ea79SLois Curfman McInnes 
13768965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
13778965ea79SLois Curfman McInnes     nz = procsnz[0];
13780752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
13798965ea79SLois Curfman McInnes 
13808965ea79SLois Curfman McInnes     /* insert into matrix */
13818965ea79SLois Curfman McInnes     jj      = rstart;
13828965ea79SLois Curfman McInnes     smycols = mycols;
13838965ea79SLois Curfman McInnes     svals   = vals;
13848965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
13858965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13868965ea79SLois Curfman McInnes       smycols += ourlens[i];
13878965ea79SLois Curfman McInnes       svals   += ourlens[i];
13888965ea79SLois Curfman McInnes       jj++;
13898965ea79SLois Curfman McInnes     }
13908965ea79SLois Curfman McInnes 
13918965ea79SLois Curfman McInnes     /* read in other processors and ship out */
13928965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
13938965ea79SLois Curfman McInnes       nz   = procsnz[i];
13940752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1395ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
13968965ea79SLois Curfman McInnes     }
1397606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
13983a40ed3dSBarry Smith   } else {
13998965ea79SLois Curfman McInnes     /* receive numeric values */
14000452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) );CHKPTRQ(vals);
14018965ea79SLois Curfman McInnes 
14028965ea79SLois Curfman McInnes     /* receive message of values*/
1403ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1404ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
1405a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
14068965ea79SLois Curfman McInnes 
14078965ea79SLois Curfman McInnes     /* insert into matrix */
14088965ea79SLois Curfman McInnes     jj      = rstart;
14098965ea79SLois Curfman McInnes     smycols = mycols;
14108965ea79SLois Curfman McInnes     svals   = vals;
14118965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
14128965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14138965ea79SLois Curfman McInnes       smycols += ourlens[i];
14148965ea79SLois Curfman McInnes       svals   += ourlens[i];
14158965ea79SLois Curfman McInnes       jj++;
14168965ea79SLois Curfman McInnes     }
14178965ea79SLois Curfman McInnes   }
1418606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1419606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1420606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1421606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
14228965ea79SLois Curfman McInnes 
14236d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14246d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14253a40ed3dSBarry Smith   PetscFunctionReturn(0);
14268965ea79SLois Curfman McInnes }
142790ace30eSBarry Smith 
142890ace30eSBarry Smith 
142990ace30eSBarry Smith 
143090ace30eSBarry Smith 
143190ace30eSBarry Smith 
1432