xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision cfce73b9bcda79ef8784c5e131f8c3eca318a602)
1*cfce73b9SSatish Balay /*$Id: mpidense.c,v 1.150 2001/01/20 03:34:41 bsmith Exp balay $*/
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__
12b0a32e0cSBarry Smith #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;
16*cfce73b9SSatish Balay   int          m = A->m,rstart = mdn->rstart,ierr;
170de54da6SSatish Balay   Scalar       *array;
180de54da6SSatish Balay   MPI_Comm     comm;
190de54da6SSatish Balay 
200de54da6SSatish Balay   PetscFunctionBegin;
21273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"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 = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
270de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
280de54da6SSatish Balay   ierr = MatCreateSeqDense(comm,m,m,array+m*rstart,B);CHKERRQ(ierr);
290de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
300de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
310de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
320de54da6SSatish Balay 
330de54da6SSatish Balay   *iscopy = PETSC_TRUE;
340de54da6SSatish Balay   PetscFunctionReturn(0);
350de54da6SSatish Balay }
360de54da6SSatish Balay EXTERN_C_END
370de54da6SSatish Balay 
38ca44d042SBarry Smith EXTERN int MatSetUpMultiply_MPIDense(Mat);
397ef1d9bdSSatish Balay 
405615d1e5SSatish Balay #undef __FUNC__
41b0a32e0cSBarry Smith #define __FUNC__ "MatSetValues_MPIDense"
428f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
438965ea79SLois Curfman McInnes {
4439b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
4539b7565bSBarry Smith   int          ierr,i,j,rstart = A->rstart,rend = A->rend,row;
46273d9f13SBarry Smith   PetscTruth   roworiented = A->roworiented;
478965ea79SLois Curfman McInnes 
483a40ed3dSBarry Smith   PetscFunctionBegin;
498965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
505ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
51273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
528965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
538965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
5439b7565bSBarry Smith       if (roworiented) {
5539b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
563a40ed3dSBarry Smith       } else {
578965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
585ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
59273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
6039b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
6139b7565bSBarry Smith         }
628965ea79SLois Curfman McInnes       }
633a40ed3dSBarry Smith     } else {
643782ba37SSatish Balay       if (!A->donotstash) {
6539b7565bSBarry Smith         if (roworiented) {
668798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
67d36fbae8SSatish Balay         } else {
688798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
6939b7565bSBarry Smith         }
70b49de8d1SLois Curfman McInnes       }
71b49de8d1SLois Curfman McInnes     }
723782ba37SSatish Balay   }
733a40ed3dSBarry Smith   PetscFunctionReturn(0);
74b49de8d1SLois Curfman McInnes }
75b49de8d1SLois Curfman McInnes 
765615d1e5SSatish Balay #undef __FUNC__
77b0a32e0cSBarry Smith #define __FUNC__ "MatGetValues_MPIDense"
788f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
79b49de8d1SLois Curfman McInnes {
80b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
81b49de8d1SLois Curfman McInnes   int          ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row;
82b49de8d1SLois Curfman McInnes 
833a40ed3dSBarry Smith   PetscFunctionBegin;
84b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
8529bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
86273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
87b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
88b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
89b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
9029bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
91273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
9229bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
93a8c6a408SBarry Smith         }
94b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
95b49de8d1SLois Curfman McInnes       }
96a8c6a408SBarry Smith     } else {
9729bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
988965ea79SLois Curfman McInnes     }
998965ea79SLois Curfman McInnes   }
1003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1018965ea79SLois Curfman McInnes }
1028965ea79SLois Curfman McInnes 
1035615d1e5SSatish Balay #undef __FUNC__
104b0a32e0cSBarry Smith #define __FUNC__ "MatGetArray_MPIDense"
1058f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array)
106ff14e315SSatish Balay {
107ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
108ff14e315SSatish Balay   int          ierr;
109ff14e315SSatish Balay 
1103a40ed3dSBarry Smith   PetscFunctionBegin;
111ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1123a40ed3dSBarry Smith   PetscFunctionReturn(0);
113ff14e315SSatish Balay }
114ff14e315SSatish Balay 
1155615d1e5SSatish Balay #undef __FUNC__
116b0a32e0cSBarry Smith #define __FUNC__ "MatGetSubMatrix_MPIDense"
117ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
118ca3fa75bSLois Curfman McInnes {
119ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
120ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
121*cfce73b9SSatish Balay   int          i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
122ca3fa75bSLois Curfman McInnes   Scalar       *av,*bv,*v = lmat->v;
123ca3fa75bSLois Curfman McInnes   Mat          newmat;
124ca3fa75bSLois Curfman McInnes 
125ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
126ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
127ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
128b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
129b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
130ca3fa75bSLois Curfman McInnes 
131ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1327eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1337eba5e9cSLois Curfman McInnes      original matrix! */
134ca3fa75bSLois Curfman McInnes 
135ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1367eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
137ca3fa75bSLois Curfman McInnes 
138ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
139ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
14029bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
1417eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
142ca3fa75bSLois Curfman McInnes     newmat = *B;
143ca3fa75bSLois Curfman McInnes   } else {
144ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
14532828cfdSBarry Smith     ierr = MatCreateMPIDense(A->comm,nrows,cs,PETSC_DECIDE,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr);
146ca3fa75bSLois Curfman McInnes   }
147ca3fa75bSLois Curfman McInnes 
148ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
149ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
150ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
151ca3fa75bSLois Curfman McInnes 
152ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
153ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
154ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
1557eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
156ca3fa75bSLois Curfman McInnes     }
157ca3fa75bSLois Curfman McInnes   }
158ca3fa75bSLois Curfman McInnes 
159ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
160ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162ca3fa75bSLois Curfman McInnes 
163ca3fa75bSLois Curfman McInnes   /* Free work space */
164ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
165ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
166ca3fa75bSLois Curfman McInnes   *B = newmat;
167ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
168ca3fa75bSLois Curfman McInnes }
169ca3fa75bSLois Curfman McInnes 
170ca3fa75bSLois Curfman McInnes #undef __FUNC__
171b0a32e0cSBarry Smith #define __FUNC__ "MatRestoreArray_MPIDense"
1728f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array)
173ff14e315SSatish Balay {
1743a40ed3dSBarry Smith   PetscFunctionBegin;
1753a40ed3dSBarry Smith   PetscFunctionReturn(0);
176ff14e315SSatish Balay }
177ff14e315SSatish Balay 
1785615d1e5SSatish Balay #undef __FUNC__
179b0a32e0cSBarry Smith #define __FUNC__ "MatAssemblyBegin_MPIDense"
1808f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
1818965ea79SLois Curfman McInnes {
18239ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
1838965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
184d36fbae8SSatish Balay   int          ierr,nstash,reallocs;
1858965ea79SLois Curfman McInnes   InsertMode   addv;
1868965ea79SLois Curfman McInnes 
1873a40ed3dSBarry Smith   PetscFunctionBegin;
1888965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
189ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
1907056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
19129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
1928965ea79SLois Curfman McInnes   }
193e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
1948965ea79SLois Curfman McInnes 
1958798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
1968798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
197b0a32e0cSBarry Smith   PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1998965ea79SLois Curfman McInnes }
2008965ea79SLois Curfman McInnes 
2015615d1e5SSatish Balay #undef __FUNC__
202b0a32e0cSBarry Smith #define __FUNC__ "MatAssemblyEnd_MPIDense"
2038f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2048965ea79SLois Curfman McInnes {
20539ddd567SLois Curfman McInnes   Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
2067ef1d9bdSSatish Balay   int          i,n,ierr,*row,*col,flg,j,rstart,ncols;
2077ef1d9bdSSatish Balay   Scalar       *val;
208e0fa3b82SLois Curfman McInnes   InsertMode   addv=mat->insertmode;
2098965ea79SLois Curfman McInnes 
2103a40ed3dSBarry Smith   PetscFunctionBegin;
2118965ea79SLois Curfman McInnes   /*  wait on receives */
2127ef1d9bdSSatish Balay   while (1) {
2138798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2147ef1d9bdSSatish Balay     if (!flg) break;
2158965ea79SLois Curfman McInnes 
2167ef1d9bdSSatish Balay     for (i=0; i<n;) {
2177ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2187ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2197ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2207ef1d9bdSSatish Balay       else       ncols = n-i;
2217ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2227ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2237ef1d9bdSSatish Balay       i = j;
2248965ea79SLois Curfman McInnes     }
2257ef1d9bdSSatish Balay   }
2268798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2278965ea79SLois Curfman McInnes 
22839ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
22939ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2308965ea79SLois Curfman McInnes 
2316d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
23239ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2338965ea79SLois Curfman McInnes   }
2343a40ed3dSBarry Smith   PetscFunctionReturn(0);
2358965ea79SLois Curfman McInnes }
2368965ea79SLois Curfman McInnes 
2375615d1e5SSatish Balay #undef __FUNC__
238b0a32e0cSBarry Smith #define __FUNC__ "MatZeroEntries_MPIDense"
2398f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A)
2408965ea79SLois Curfman McInnes {
2413a40ed3dSBarry Smith   int          ierr;
24239ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense*)A->data;
2433a40ed3dSBarry Smith 
2443a40ed3dSBarry Smith   PetscFunctionBegin;
2453a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2463a40ed3dSBarry Smith   PetscFunctionReturn(0);
2478965ea79SLois Curfman McInnes }
2488965ea79SLois Curfman McInnes 
2495615d1e5SSatish Balay #undef __FUNC__
250b0a32e0cSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIDense"
2518f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs)
2524e220ebcSLois Curfman McInnes {
2533a40ed3dSBarry Smith   PetscFunctionBegin;
2544e220ebcSLois Curfman McInnes   *bs = 1;
2553a40ed3dSBarry Smith   PetscFunctionReturn(0);
2564e220ebcSLois Curfman McInnes }
2574e220ebcSLois Curfman McInnes 
2588965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2598965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2608965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2613501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2628965ea79SLois Curfman McInnes    routine.
2638965ea79SLois Curfman McInnes */
2645615d1e5SSatish Balay #undef __FUNC__
265b0a32e0cSBarry Smith #define __FUNC__ "MatZeroRows_MPIDense"
2668f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag)
2678965ea79SLois Curfman McInnes {
26839ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
2698965ea79SLois Curfman McInnes   int            i,ierr,N,*rows,*owners = l->rowners,size = l->size;
27035d8aa7fSBarry Smith   int            *procs,*nprocs,j,idx,nsends,*work;
2718965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2728965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2738965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2748965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2758965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2768965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2778965ea79SLois Curfman McInnes   IS             istmp;
27835d8aa7fSBarry Smith   PetscTruth     found;
2798965ea79SLois Curfman McInnes 
2803a40ed3dSBarry Smith   PetscFunctionBegin;
281b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
2828965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
2838965ea79SLois Curfman McInnes 
2848965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
285b0a32e0cSBarry Smith   ierr  = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr);
286549d3d68SSatish Balay   ierr  = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
287549d3d68SSatish Balay   procs = nprocs + size;
288b0a32e0cSBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/
2898965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
2908965ea79SLois Curfman McInnes     idx = rows[i];
29135d8aa7fSBarry Smith     found = PETSC_FALSE;
2928965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
2938965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
29435d8aa7fSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; found = PETSC_TRUE; break;
2958965ea79SLois Curfman McInnes       }
2968965ea79SLois Curfman McInnes     }
29729bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
2988965ea79SLois Curfman McInnes   }
2998965ea79SLois Curfman McInnes   nsends = 0;  for (i=0; i<size; i++) { nsends += procs[i];}
3008965ea79SLois Curfman McInnes 
3018965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
302b0a32e0cSBarry Smith   ierr   = PetscMalloc(2*size*sizeof(int),&work);CHKERRQ(ierr);
3036831982aSBarry Smith   ierr   = MPI_Allreduce(nprocs,work,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr);
3048965ea79SLois Curfman McInnes   nmax   = work[rank];
3056831982aSBarry Smith   nrecvs = work[size+rank];
306606d414cSSatish Balay   ierr   = PetscFree(work);CHKERRQ(ierr);
3078965ea79SLois Curfman McInnes 
3088965ea79SLois Curfman McInnes   /* post receives:   */
309b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr);
310b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3118965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
312ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3138965ea79SLois Curfman McInnes   }
3148965ea79SLois Curfman McInnes 
3158965ea79SLois Curfman McInnes   /* do sends:
3168965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3178965ea79SLois Curfman McInnes          the ith processor
3188965ea79SLois Curfman McInnes   */
319b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr);
320b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
321b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr);
3228965ea79SLois Curfman McInnes   starts[0]  = 0;
3238965ea79SLois Curfman McInnes   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];}
3248965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3258965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3268965ea79SLois Curfman McInnes   }
3278965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3288965ea79SLois Curfman McInnes 
3298965ea79SLois Curfman McInnes   starts[0] = 0;
3308965ea79SLois Curfman McInnes   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[i-1];}
3318965ea79SLois Curfman McInnes   count = 0;
3328965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
3338965ea79SLois Curfman McInnes     if (procs[i]) {
334ca161407SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3358965ea79SLois Curfman McInnes     }
3368965ea79SLois Curfman McInnes   }
337606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3388965ea79SLois Curfman McInnes 
3398965ea79SLois Curfman McInnes   base = owners[rank];
3408965ea79SLois Curfman McInnes 
3418965ea79SLois Curfman McInnes   /*  wait on receives */
342b0a32e0cSBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr);
3438965ea79SLois Curfman McInnes   source = lens + nrecvs;
3448965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3458965ea79SLois Curfman McInnes   while (count) {
346ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3478965ea79SLois Curfman McInnes     /* unpack receives into our local space */
348ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3498965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3508965ea79SLois Curfman McInnes     lens[imdex]  = n;
3518965ea79SLois Curfman McInnes     slen += n;
3528965ea79SLois Curfman McInnes     count--;
3538965ea79SLois Curfman McInnes   }
354606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3558965ea79SLois Curfman McInnes 
3568965ea79SLois Curfman McInnes   /* move the data into the send scatter */
357b0a32e0cSBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr);
3588965ea79SLois Curfman McInnes   count = 0;
3598965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
3608965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3618965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
3628965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3638965ea79SLois Curfman McInnes     }
3648965ea79SLois Curfman McInnes   }
365606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
366606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
367606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
368606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3698965ea79SLois Curfman McInnes 
3708965ea79SLois Curfman McInnes   /* actually zap the local rows */
371029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
372b0a32e0cSBarry Smith   PetscLogObjectParent(A,istmp);
373606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3748965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3758965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3768965ea79SLois Curfman McInnes 
3778965ea79SLois Curfman McInnes   /* wait on sends */
3788965ea79SLois Curfman McInnes   if (nsends) {
379b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
380ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
381606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
3828965ea79SLois Curfman McInnes   }
383606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
384606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
3858965ea79SLois Curfman McInnes 
3863a40ed3dSBarry Smith   PetscFunctionReturn(0);
3878965ea79SLois Curfman McInnes }
3888965ea79SLois Curfman McInnes 
3895615d1e5SSatish Balay #undef __FUNC__
390b0a32e0cSBarry Smith #define __FUNC__ "MatMult_MPIDense"
3918f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3928965ea79SLois Curfman McInnes {
39339ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
3948965ea79SLois Curfman McInnes   int          ierr;
395c456f294SBarry Smith 
3963a40ed3dSBarry Smith   PetscFunctionBegin;
39743a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39843a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39944cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4003a40ed3dSBarry Smith   PetscFunctionReturn(0);
4018965ea79SLois Curfman McInnes }
4028965ea79SLois Curfman McInnes 
4035615d1e5SSatish Balay #undef __FUNC__
404b0a32e0cSBarry Smith #define __FUNC__ "MatMultAdd_MPIDense"
4058f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4068965ea79SLois Curfman McInnes {
40739ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4088965ea79SLois Curfman McInnes   int          ierr;
409c456f294SBarry Smith 
4103a40ed3dSBarry Smith   PetscFunctionBegin;
41143a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41243a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41344cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4143a40ed3dSBarry Smith   PetscFunctionReturn(0);
4158965ea79SLois Curfman McInnes }
4168965ea79SLois Curfman McInnes 
4175615d1e5SSatish Balay #undef __FUNC__
418b0a32e0cSBarry Smith #define __FUNC__ "MatMultTranspose_MPIDense"
4197c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
420096963f5SLois Curfman McInnes {
421096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
422096963f5SLois Curfman McInnes   int          ierr;
4233501a2bdSLois Curfman McInnes   Scalar       zero = 0.0;
424096963f5SLois Curfman McInnes 
4253a40ed3dSBarry Smith   PetscFunctionBegin;
4263501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
4277c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
428537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
429537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4303a40ed3dSBarry Smith   PetscFunctionReturn(0);
431096963f5SLois Curfman McInnes }
432096963f5SLois Curfman McInnes 
4335615d1e5SSatish Balay #undef __FUNC__
434b0a32e0cSBarry Smith #define __FUNC__ "MatMultTransposeAdd_MPIDense"
4357c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
436096963f5SLois Curfman McInnes {
437096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
438096963f5SLois Curfman McInnes   int          ierr;
439096963f5SLois Curfman McInnes 
4403a40ed3dSBarry Smith   PetscFunctionBegin;
4413501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4427c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
443537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
444537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4453a40ed3dSBarry Smith   PetscFunctionReturn(0);
446096963f5SLois Curfman McInnes }
447096963f5SLois Curfman McInnes 
4485615d1e5SSatish Balay #undef __FUNC__
449b0a32e0cSBarry Smith #define __FUNC__ "MatGetDiagonal_MPIDense"
4508f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v)
4518965ea79SLois Curfman McInnes {
45239ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
453096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
454273d9f13SBarry Smith   int          ierr,len,i,n,m = A->m,radd;
45544cd7ae7SLois Curfman McInnes   Scalar       *x,zero = 0.0;
456ed3cc1f0SBarry Smith 
4573a40ed3dSBarry Smith   PetscFunctionBegin;
458273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
459096963f5SLois Curfman McInnes   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
460096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
461273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
462273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
4637ddc982cSLois Curfman McInnes   radd = a->rstart*m;
46444cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
465096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
466096963f5SLois Curfman McInnes   }
4679a8c540fSSatish Balay   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
4683a40ed3dSBarry Smith   PetscFunctionReturn(0);
4698965ea79SLois Curfman McInnes }
4708965ea79SLois Curfman McInnes 
4715615d1e5SSatish Balay #undef __FUNC__
472b0a32e0cSBarry Smith #define __FUNC__ "MatDestroy_MPIDense"
473e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat)
4748965ea79SLois Curfman McInnes {
4753501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4768965ea79SLois Curfman McInnes   int          ierr;
477ed3cc1f0SBarry Smith 
4783a40ed3dSBarry Smith   PetscFunctionBegin;
47994d884c6SBarry Smith 
480aa482453SBarry Smith #if defined(PETSC_USE_LOG)
481b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N);
4828965ea79SLois Curfman McInnes #endif
4838798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
484606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
4853501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
4863501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4873501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
488622d7880SLois Curfman McInnes   if (mdn->factor) {
489606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
490606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
491606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
492606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
493622d7880SLois Curfman McInnes   }
494606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
4953a40ed3dSBarry Smith   PetscFunctionReturn(0);
4968965ea79SLois Curfman McInnes }
49739ddd567SLois Curfman McInnes 
4985615d1e5SSatish Balay #undef __FUNC__
499b0a32e0cSBarry Smith #define __FUNC__ "MatView_MPIDense_Binary"
500b0a32e0cSBarry Smith static int MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5018965ea79SLois Curfman McInnes {
50239ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
5038965ea79SLois Curfman McInnes   int          ierr;
5047056b6fcSBarry Smith 
5053a40ed3dSBarry Smith   PetscFunctionBegin;
50639ddd567SLois Curfman McInnes   if (mdn->size == 1) {
50739ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5088965ea79SLois Curfman McInnes   }
50929bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5103a40ed3dSBarry Smith   PetscFunctionReturn(0);
5118965ea79SLois Curfman McInnes }
5128965ea79SLois Curfman McInnes 
5135615d1e5SSatish Balay #undef __FUNC__
514b0a32e0cSBarry Smith #define __FUNC__ "MatView_MPIDense_ASCIIorDraworSocket"
515b0a32e0cSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5168965ea79SLois Curfman McInnes {
51739ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
518fb9695e5SSatish Balay   int               ierr,size = mdn->size,rank = mdn->rank;
519b0a32e0cSBarry Smith   PetscViewerType   vtype;
520f1af5d2fSBarry Smith   PetscTruth        isascii,isdraw;
521b0a32e0cSBarry Smith   PetscViewer       sviewer;
522f3ef73ceSBarry Smith   PetscViewerFormat format;
5238965ea79SLois Curfman McInnes 
5243a40ed3dSBarry Smith   PetscFunctionBegin;
525b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
526fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
527f1af5d2fSBarry Smith   if (isascii) {
528b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
529b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
530fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_INFO_LONG) {
5314e220ebcSLois Curfman McInnes       MatInfo info;
532888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
533b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m,
5346831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
535b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5363501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5373a40ed3dSBarry Smith       PetscFunctionReturn(0);
538fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5393a40ed3dSBarry Smith       PetscFunctionReturn(0);
5408965ea79SLois Curfman McInnes     }
541f1af5d2fSBarry Smith   } else if (isdraw) {
542b0a32e0cSBarry Smith     PetscDraw       draw;
543f1af5d2fSBarry Smith     PetscTruth isnull;
544f1af5d2fSBarry Smith 
545b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
546b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
547f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
548f1af5d2fSBarry Smith   }
54977ed5343SBarry Smith 
5508965ea79SLois Curfman McInnes   if (size == 1) {
55139ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5523a40ed3dSBarry Smith   } else {
5538965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5548965ea79SLois Curfman McInnes     Mat          A;
555273d9f13SBarry Smith     int          M = mat->M,N = mat->N,m,row,i,nz,*cols;
55639ddd567SLois Curfman McInnes     Scalar       *vals;
5578965ea79SLois Curfman McInnes 
5588965ea79SLois Curfman McInnes     if (!rank) {
559f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5603a40ed3dSBarry Smith     } else {
561f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5628965ea79SLois Curfman McInnes     }
563b0a32e0cSBarry Smith     PetscLogObjectParent(mat,A);
5648965ea79SLois Curfman McInnes 
56539ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
56639ddd567SLois Curfman McInnes        but it's quick for now */
567273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
5688965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
56939ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57039ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
57139ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57239ddd567SLois Curfman McInnes       row++;
5738965ea79SLois Curfman McInnes     }
5748965ea79SLois Curfman McInnes 
5756d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5766d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
577b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
578b9b97703SBarry Smith     if (!rank) {
5796831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5808965ea79SLois Curfman McInnes     }
581b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
582b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5838965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
5848965ea79SLois Curfman McInnes   }
5853a40ed3dSBarry Smith   PetscFunctionReturn(0);
5868965ea79SLois Curfman McInnes }
5878965ea79SLois Curfman McInnes 
5885615d1e5SSatish Balay #undef __FUNC__
589b0a32e0cSBarry Smith #define __FUNC__ "MatView_MPIDense"
590b0a32e0cSBarry Smith int MatView_MPIDense(Mat mat,PetscViewer viewer)
5918965ea79SLois Curfman McInnes {
59239ddd567SLois Curfman McInnes   int        ierr;
593f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
5948965ea79SLois Curfman McInnes 
595433994e6SBarry Smith   PetscFunctionBegin;
5960f5bd95cSBarry Smith 
597b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
598fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
599b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
600fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6010f5bd95cSBarry Smith 
602f1af5d2fSBarry Smith   if (isascii || issocket || isdraw) {
603f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6040f5bd95cSBarry Smith   } else if (isbinary) {
6053a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6065cd90555SBarry Smith   } else {
60729bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6088965ea79SLois Curfman McInnes   }
6093a40ed3dSBarry Smith   PetscFunctionReturn(0);
6108965ea79SLois Curfman McInnes }
6118965ea79SLois Curfman McInnes 
6125615d1e5SSatish Balay #undef __FUNC__
613b0a32e0cSBarry Smith #define __FUNC__ "MatGetInfo_MPIDense"
6148f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6158965ea79SLois Curfman McInnes {
6163501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6173501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6184e220ebcSLois Curfman McInnes   int          ierr;
619329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6208965ea79SLois Curfman McInnes 
6213a40ed3dSBarry Smith   PetscFunctionBegin;
622273d9f13SBarry Smith   info->rows_global    = (double)A->M;
623273d9f13SBarry Smith   info->columns_global = (double)A->N;
624273d9f13SBarry Smith   info->rows_local     = (double)A->m;
625273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6264e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6274e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6284e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6294e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6308965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6314e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6324e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6334e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6344e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6354e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6368965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
637f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
6384e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6394e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6404e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6414e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6424e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6438965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
644f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
6454e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6464e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6474e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6484e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6494e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6508965ea79SLois Curfman McInnes   }
6514e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6524e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6534e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6543a40ed3dSBarry Smith   PetscFunctionReturn(0);
6558965ea79SLois Curfman McInnes }
6568965ea79SLois Curfman McInnes 
6575615d1e5SSatish Balay #undef __FUNC__
658b0a32e0cSBarry Smith #define __FUNC__ "MatSetOption_MPIDense"
6598f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6608965ea79SLois Curfman McInnes {
66139ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
662273d9f13SBarry Smith   int          ierr;
6638965ea79SLois Curfman McInnes 
6643a40ed3dSBarry Smith   PetscFunctionBegin;
6656d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6666d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
6674787f768SSatish Balay       op == MAT_NEW_NONZERO_LOCATION_ERR ||
6684787f768SSatish Balay       op == MAT_NEW_NONZERO_ALLOCATION_ERR ||
669219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
670219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED) {
671273d9f13SBarry Smith         ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
672b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
673273d9f13SBarry Smith         a->roworiented = PETSC_TRUE;
674273d9f13SBarry Smith         ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
675b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
676219d9a1aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
6776d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
6786d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
679b51ba29fSSatish Balay              op == MAT_YES_NEW_DIAGONALS ||
680b51ba29fSSatish Balay              op == MAT_USE_HASH_TABLE) {
681b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
6823a40ed3dSBarry Smith   } else if (op == MAT_COLUMN_ORIENTED) {
683273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
684273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
6853782ba37SSatish Balay   } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) {
686273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
6873a40ed3dSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS) {
68829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
6893a40ed3dSBarry Smith   } else {
69029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
6913a40ed3dSBarry Smith   }
6923a40ed3dSBarry Smith   PetscFunctionReturn(0);
6938965ea79SLois Curfman McInnes }
6948965ea79SLois Curfman McInnes 
6955615d1e5SSatish Balay #undef __FUNC__
696b0a32e0cSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIDense"
6978f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
6988965ea79SLois Curfman McInnes {
6993501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7003a40ed3dSBarry Smith 
7013a40ed3dSBarry Smith   PetscFunctionBegin;
7024c49b128SBarry Smith   if (m) *m = mat->rstart;
7034c49b128SBarry Smith   if (n) *n = mat->rend;
7043a40ed3dSBarry Smith   PetscFunctionReturn(0);
7058965ea79SLois Curfman McInnes }
7068965ea79SLois Curfman McInnes 
7075615d1e5SSatish Balay #undef __FUNC__
708b0a32e0cSBarry Smith #define __FUNC__ "MatGetRow_MPIDense"
7098f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
7108965ea79SLois Curfman McInnes {
7113501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7123a40ed3dSBarry Smith   int          lrow,rstart = mat->rstart,rend = mat->rend,ierr;
7138965ea79SLois Curfman McInnes 
7143a40ed3dSBarry Smith   PetscFunctionBegin;
71529bbc08cSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
7168965ea79SLois Curfman McInnes   lrow = row - rstart;
7173a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7183a40ed3dSBarry Smith   PetscFunctionReturn(0);
7198965ea79SLois Curfman McInnes }
7208965ea79SLois Curfman McInnes 
7215615d1e5SSatish Balay #undef __FUNC__
722b0a32e0cSBarry Smith #define __FUNC__ "MatRestoreRow_MPIDense"
7238f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
7248965ea79SLois Curfman McInnes {
725606d414cSSatish Balay   int ierr;
726606d414cSSatish Balay 
7273a40ed3dSBarry Smith   PetscFunctionBegin;
728606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
729606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7303a40ed3dSBarry Smith   PetscFunctionReturn(0);
7318965ea79SLois Curfman McInnes }
7328965ea79SLois Curfman McInnes 
7335615d1e5SSatish Balay #undef __FUNC__
734b0a32e0cSBarry Smith #define __FUNC__ "MatDiagonalScale_MPIDense"
7355b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7365b2fa520SLois Curfman McInnes {
7375b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7385b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7395b2fa520SLois Curfman McInnes   Scalar       *l,*r,x,*v;
740273d9f13SBarry Smith   int          ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7415b2fa520SLois Curfman McInnes 
7425b2fa520SLois Curfman McInnes   PetscFunctionBegin;
74372d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7445b2fa520SLois Curfman McInnes   if (ll) {
74572d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
74629bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
7475b2fa520SLois Curfman McInnes     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7485b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7495b2fa520SLois Curfman McInnes       x = l[i];
7505b2fa520SLois Curfman McInnes       v = mat->v + i;
7515b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7525b2fa520SLois Curfman McInnes     }
7535b2fa520SLois Curfman McInnes     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
754b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7555b2fa520SLois Curfman McInnes   }
7565b2fa520SLois Curfman McInnes   if (rr) {
75772d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
75829bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7595b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7605b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7615b2fa520SLois Curfman McInnes     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
7625b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7635b2fa520SLois Curfman McInnes       x = r[i];
7645b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7655b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7665b2fa520SLois Curfman McInnes     }
76772d926a5SLois Curfman McInnes     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
768b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7695b2fa520SLois Curfman McInnes   }
7705b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7715b2fa520SLois Curfman McInnes }
7725b2fa520SLois Curfman McInnes 
7735b2fa520SLois Curfman McInnes #undef __FUNC__
774b0a32e0cSBarry Smith #define __FUNC__ "MatNorm_MPIDense"
775329f5518SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *norm)
776096963f5SLois Curfman McInnes {
7773501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7783501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7793501a2bdSLois Curfman McInnes   int          ierr,i,j;
780329f5518SBarry Smith   PetscReal    sum = 0.0;
7813501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
7823501a2bdSLois Curfman McInnes 
7833a40ed3dSBarry Smith   PetscFunctionBegin;
7843501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
7853501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm);CHKERRQ(ierr);
7863501a2bdSLois Curfman McInnes   } else {
7873501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
788273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
789aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
790329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
7913501a2bdSLois Curfman McInnes #else
7923501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7933501a2bdSLois Curfman McInnes #endif
7943501a2bdSLois Curfman McInnes       }
795ca161407SBarry Smith       ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
7963501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
797b0a32e0cSBarry Smith       PetscLogFlops(2*mdn->A->n*mdn->A->m);
7983a40ed3dSBarry Smith     } else if (type == NORM_1) {
799329f5518SBarry Smith       PetscReal *tmp,*tmp2;
800b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
801273d9f13SBarry Smith       tmp2 = tmp + A->N;
802273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
803096963f5SLois Curfman McInnes       *norm = 0.0;
8043501a2bdSLois Curfman McInnes       v = mat->v;
805273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
806273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
80767e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8083501a2bdSLois Curfman McInnes         }
8093501a2bdSLois Curfman McInnes       }
810273d9f13SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
811273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
8123501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
8133501a2bdSLois Curfman McInnes       }
814606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
815b0a32e0cSBarry Smith       PetscLogFlops(A->n*A->m);
8163a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
817329f5518SBarry Smith       PetscReal ntemp;
8183501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
819ca161407SBarry Smith       ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
8203a40ed3dSBarry Smith     } else {
82129bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8223501a2bdSLois Curfman McInnes     }
8233501a2bdSLois Curfman McInnes   }
8243a40ed3dSBarry Smith   PetscFunctionReturn(0);
8253501a2bdSLois Curfman McInnes }
8263501a2bdSLois Curfman McInnes 
8275615d1e5SSatish Balay #undef __FUNC__
828b0a32e0cSBarry Smith #define __FUNC__ "MatTranspose_MPIDense"
8298f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8303501a2bdSLois Curfman McInnes {
8313501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8323501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8333501a2bdSLois Curfman McInnes   Mat          B;
834273d9f13SBarry Smith   int          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8353501a2bdSLois Curfman McInnes   int          j,i,ierr;
8363501a2bdSLois Curfman McInnes   Scalar       *v;
8373501a2bdSLois Curfman McInnes 
8383a40ed3dSBarry Smith   PetscFunctionBegin;
8397c922b88SBarry Smith   if (!matout && M != N) {
84029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8417056b6fcSBarry Smith   }
8427056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8433501a2bdSLois Curfman McInnes 
844273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
845b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr);
8463501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8473501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8483501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8493501a2bdSLois Curfman McInnes     v   += m;
8503501a2bdSLois Curfman McInnes   }
851606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8526d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8536d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8547c922b88SBarry Smith   if (matout) {
8553501a2bdSLois Curfman McInnes     *matout = B;
8563501a2bdSLois Curfman McInnes   } else {
857273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8583501a2bdSLois Curfman McInnes   }
8593a40ed3dSBarry Smith   PetscFunctionReturn(0);
860096963f5SLois Curfman McInnes }
861096963f5SLois Curfman McInnes 
862d9eff348SSatish Balay #include "petscblaslapack.h"
8635615d1e5SSatish Balay #undef __FUNC__
864b0a32e0cSBarry Smith #define __FUNC__ "MatScale_MPIDense"
8658f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA)
86644cd7ae7SLois Curfman McInnes {
86744cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
86844cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
86944cd7ae7SLois Curfman McInnes   int          one = 1,nz;
87044cd7ae7SLois Curfman McInnes 
8713a40ed3dSBarry Smith   PetscFunctionBegin;
872273d9f13SBarry Smith   nz = inA->m*inA->N;
87344cd7ae7SLois Curfman McInnes   BLscal_(&nz,alpha,a->v,&one);
874b0a32e0cSBarry Smith   PetscLogFlops(nz);
8753a40ed3dSBarry Smith   PetscFunctionReturn(0);
87644cd7ae7SLois Curfman McInnes }
87744cd7ae7SLois Curfman McInnes 
8785609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
879ca44d042SBarry Smith EXTERN int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **);
8808965ea79SLois Curfman McInnes 
881273d9f13SBarry Smith #undef __FUNC__
882b0a32e0cSBarry Smith #define __FUNC__ "MatSetUpPreallocation_MPIDense"
883273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A)
884273d9f13SBarry Smith {
885273d9f13SBarry Smith   int        ierr;
886273d9f13SBarry Smith 
887273d9f13SBarry Smith   PetscFunctionBegin;
888273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
889273d9f13SBarry Smith   PetscFunctionReturn(0);
890273d9f13SBarry Smith }
891273d9f13SBarry Smith 
8928965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
89309dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
89409dc0095SBarry Smith        MatGetRow_MPIDense,
89509dc0095SBarry Smith        MatRestoreRow_MPIDense,
89609dc0095SBarry Smith        MatMult_MPIDense,
89709dc0095SBarry Smith        MatMultAdd_MPIDense,
8987c922b88SBarry Smith        MatMultTranspose_MPIDense,
8997c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9008965ea79SLois Curfman McInnes        0,
90109dc0095SBarry Smith        0,
90209dc0095SBarry Smith        0,
90309dc0095SBarry Smith        0,
90409dc0095SBarry Smith        0,
90509dc0095SBarry Smith        0,
90609dc0095SBarry Smith        0,
90709dc0095SBarry Smith        MatTranspose_MPIDense,
90809dc0095SBarry Smith        MatGetInfo_MPIDense,0,
90909dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9105b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
91109dc0095SBarry Smith        MatNorm_MPIDense,
91209dc0095SBarry Smith        MatAssemblyBegin_MPIDense,
91309dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
91409dc0095SBarry Smith        0,
91509dc0095SBarry Smith        MatSetOption_MPIDense,
91609dc0095SBarry Smith        MatZeroEntries_MPIDense,
91709dc0095SBarry Smith        MatZeroRows_MPIDense,
91809dc0095SBarry Smith        0,
91909dc0095SBarry Smith        0,
92009dc0095SBarry Smith        0,
92109dc0095SBarry Smith        0,
922273d9f13SBarry Smith        MatSetUpPreallocation_MPIDense,
923273d9f13SBarry Smith        0,
92439ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
92509dc0095SBarry Smith        0,
92609dc0095SBarry Smith        0,
92709dc0095SBarry Smith        MatGetArray_MPIDense,
92809dc0095SBarry Smith        MatRestoreArray_MPIDense,
9295609ef8eSBarry Smith        MatDuplicate_MPIDense,
93009dc0095SBarry Smith        0,
93109dc0095SBarry Smith        0,
93209dc0095SBarry Smith        0,
93309dc0095SBarry Smith        0,
93409dc0095SBarry Smith        0,
9352ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
93609dc0095SBarry Smith        0,
93709dc0095SBarry Smith        MatGetValues_MPIDense,
93809dc0095SBarry Smith        0,
93909dc0095SBarry Smith        0,
94009dc0095SBarry Smith        MatScale_MPIDense,
94109dc0095SBarry Smith        0,
94209dc0095SBarry Smith        0,
94309dc0095SBarry Smith        0,
94409dc0095SBarry Smith        MatGetBlockSize_MPIDense,
94509dc0095SBarry Smith        0,
94609dc0095SBarry Smith        0,
94709dc0095SBarry Smith        0,
94809dc0095SBarry Smith        0,
94909dc0095SBarry Smith        0,
95009dc0095SBarry Smith        0,
95109dc0095SBarry Smith        0,
95209dc0095SBarry Smith        0,
95309dc0095SBarry Smith        0,
954ca3fa75bSLois Curfman McInnes        MatGetSubMatrix_MPIDense,
955b9b97703SBarry Smith        MatDestroy_MPIDense,
956b9b97703SBarry Smith        MatView_MPIDense,
95709dc0095SBarry Smith        MatGetMaps_Petsc};
9588965ea79SLois Curfman McInnes 
959273d9f13SBarry Smith EXTERN_C_BEGIN
960273d9f13SBarry Smith #undef __FUNC__
961b0a32e0cSBarry Smith #define __FUNC__ "MatCreate_MPIDense"
962273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat)
963273d9f13SBarry Smith {
964273d9f13SBarry Smith   Mat_MPIDense *a;
965273d9f13SBarry Smith   int          ierr,i;
966273d9f13SBarry Smith 
967273d9f13SBarry Smith   PetscFunctionBegin;
968b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
969b0a32e0cSBarry Smith   mat->data         = (void*)a;
970273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
971273d9f13SBarry Smith   mat->factor       = 0;
972273d9f13SBarry Smith   mat->mapping      = 0;
973273d9f13SBarry Smith 
974273d9f13SBarry Smith   a->factor       = 0;
975273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
976273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
977273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
978273d9f13SBarry Smith 
979273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
980273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
981273d9f13SBarry Smith   a->nvec = mat->n;
982273d9f13SBarry Smith 
983273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
984273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
985273d9f13SBarry Smith   ierr = MapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
986273d9f13SBarry Smith   ierr = MapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
987273d9f13SBarry Smith 
988273d9f13SBarry Smith   /* build local table of row and column ownerships */
989b0a32e0cSBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr);
990273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
991b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
992273d9f13SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
993273d9f13SBarry Smith   a->rowners[0] = 0;
994273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
995273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
996273d9f13SBarry Smith   }
997273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
998273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
999273d9f13SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1000273d9f13SBarry Smith   a->cowners[0] = 0;
1001273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1002273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1003273d9f13SBarry Smith   }
1004273d9f13SBarry Smith 
1005273d9f13SBarry Smith   /* build cache for off array entries formed */
1006273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1007273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1008273d9f13SBarry Smith 
1009273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1010273d9f13SBarry Smith   a->lvec        = 0;
1011273d9f13SBarry Smith   a->Mvctx       = 0;
1012273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1013273d9f13SBarry Smith 
1014273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1015273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1016273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1017273d9f13SBarry Smith   PetscFunctionReturn(0);
1018273d9f13SBarry Smith }
1019273d9f13SBarry Smith EXTERN_C_END
1020273d9f13SBarry Smith 
1021273d9f13SBarry Smith #undef __FUNC__
1022b0a32e0cSBarry Smith #define __FUNC__ "MatMPIDenseSetPreallocation"
1023273d9f13SBarry Smith /*@C
1024273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1025273d9f13SBarry Smith 
1026273d9f13SBarry Smith    Not collective
1027273d9f13SBarry Smith 
1028273d9f13SBarry Smith    Input Parameters:
1029273d9f13SBarry Smith .  A - the matrix
1030273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1031273d9f13SBarry Smith    to control all matrix memory allocation.
1032273d9f13SBarry Smith 
1033273d9f13SBarry Smith    Notes:
1034273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1035273d9f13SBarry Smith    storage by columns.
1036273d9f13SBarry Smith 
1037273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1038273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1039273d9f13SBarry Smith    set data=PETSC_NULL.
1040273d9f13SBarry Smith 
1041273d9f13SBarry Smith    Level: intermediate
1042273d9f13SBarry Smith 
1043273d9f13SBarry Smith .keywords: matrix,dense, parallel
1044273d9f13SBarry Smith 
1045273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1046273d9f13SBarry Smith @*/
1047273d9f13SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,Scalar *data)
1048273d9f13SBarry Smith {
1049273d9f13SBarry Smith   Mat_MPIDense *a;
1050273d9f13SBarry Smith   int          ierr;
1051273d9f13SBarry Smith   PetscTruth   flg2;
1052273d9f13SBarry Smith 
1053273d9f13SBarry Smith   PetscFunctionBegin;
1054273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)mat,MATMPIDENSE,&flg2);CHKERRQ(ierr);
1055273d9f13SBarry Smith   if (!flg2) PetscFunctionReturn(0);
1056273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
1057273d9f13SBarry Smith   /* Note:  For now, when data is specified above, this assumes the user correctly
1058273d9f13SBarry Smith    allocates the local dense storage space.  We should add error checking. */
1059273d9f13SBarry Smith 
1060273d9f13SBarry Smith   a    = (Mat_MPIDense*)mat->data;
1061273d9f13SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mat->m,mat->N,data,&a->A);CHKERRQ(ierr);
1062b0a32e0cSBarry Smith   PetscLogObjectParent(mat,a->A);
1063273d9f13SBarry Smith   PetscFunctionReturn(0);
1064273d9f13SBarry Smith }
1065273d9f13SBarry Smith 
10665615d1e5SSatish Balay #undef __FUNC__
1067b0a32e0cSBarry Smith #define __FUNC__ "MatCreateMPIDense"
10688965ea79SLois Curfman McInnes /*@C
106939ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
10708965ea79SLois Curfman McInnes 
1071db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1072db81eaa0SLois Curfman McInnes 
10738965ea79SLois Curfman McInnes    Input Parameters:
1074db81eaa0SLois Curfman McInnes +  comm - MPI communicator
10758965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1076db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
10778965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1078db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1079db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1080dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
10818965ea79SLois Curfman McInnes 
10828965ea79SLois Curfman McInnes    Output Parameter:
1083477f1c0bSLois Curfman McInnes .  A - the matrix
10848965ea79SLois Curfman McInnes 
1085b259b22eSLois Curfman McInnes    Notes:
108639ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
108739ddd567SLois Curfman McInnes    storage by columns.
10888965ea79SLois Curfman McInnes 
108918f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
109018f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1091b4fd4287SBarry Smith    set data=PETSC_NULL.
109218f449edSLois Curfman McInnes 
10938965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
10948965ea79SLois Curfman McInnes    (possibly both).
10958965ea79SLois Curfman McInnes 
1096027ccd11SLois Curfman McInnes    Level: intermediate
1097027ccd11SLois Curfman McInnes 
109839ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
10998965ea79SLois Curfman McInnes 
110039ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
11018965ea79SLois Curfman McInnes @*/
1102477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
11038965ea79SLois Curfman McInnes {
1104273d9f13SBarry Smith   int ierr,size;
11058965ea79SLois Curfman McInnes 
11063a40ed3dSBarry Smith   PetscFunctionBegin;
1107273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1108273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1109273d9f13SBarry Smith   if (size > 1) {
1110273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1111273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1112273d9f13SBarry Smith   } else {
1113273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1114273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
11158c469469SLois Curfman McInnes   }
11163a40ed3dSBarry Smith   PetscFunctionReturn(0);
11178965ea79SLois Curfman McInnes }
11188965ea79SLois Curfman McInnes 
11195615d1e5SSatish Balay #undef __FUNC__
1120b0a32e0cSBarry Smith #define __FUNC__ "MatDuplicate_MPIDense"
11215609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11228965ea79SLois Curfman McInnes {
11238965ea79SLois Curfman McInnes   Mat          mat;
11243501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
112539ddd567SLois Curfman McInnes   int          ierr;
11268965ea79SLois Curfman McInnes 
11273a40ed3dSBarry Smith   PetscFunctionBegin;
11288965ea79SLois Curfman McInnes   *newmat       = 0;
1129273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1130273d9f13SBarry Smith   ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr);
1131b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1132b0a32e0cSBarry Smith   mat->data         = (void*)a;
1133549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
11343501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1135c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1136273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
11378965ea79SLois Curfman McInnes 
11388965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
11398965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
11408965ea79SLois Curfman McInnes   a->size         = oldmat->size;
11418965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1142e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1143b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
11443782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1145b0a32e0cSBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1146b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1147549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
11488798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
11498965ea79SLois Curfman McInnes 
1150329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
11515609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
1152b0a32e0cSBarry Smith   PetscLogObjectParent(mat,a->A);
11538965ea79SLois Curfman McInnes   *newmat = mat;
11543a40ed3dSBarry Smith   PetscFunctionReturn(0);
11558965ea79SLois Curfman McInnes }
11568965ea79SLois Curfman McInnes 
1157e090d566SSatish Balay #include "petscsys.h"
11588965ea79SLois Curfman McInnes 
11595615d1e5SSatish Balay #undef __FUNC__
1160b0a32e0cSBarry Smith #define __FUNC__ "MatLoad_MPIDense_DenseInFile"
116190ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat)
116290ace30eSBarry Smith {
116340011551SBarry Smith   int        *rowners,i,size,rank,m,ierr,nz,j;
116490ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
116590ace30eSBarry Smith   MPI_Status status;
116690ace30eSBarry Smith 
11673a40ed3dSBarry Smith   PetscFunctionBegin;
1168d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1169d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
117090ace30eSBarry Smith 
117190ace30eSBarry Smith   /* determine ownership of all rows */
117290ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
1173b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1174ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
117590ace30eSBarry Smith   rowners[0] = 0;
117690ace30eSBarry Smith   for (i=2; i<=size; i++) {
117790ace30eSBarry Smith     rowners[i] += rowners[i-1];
117890ace30eSBarry Smith   }
117990ace30eSBarry Smith 
118090ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
118190ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
118290ace30eSBarry Smith 
118390ace30eSBarry Smith   if (!rank) {
1184b0a32e0cSBarry Smith     ierr = PetscMalloc(m*N*sizeof(Scalar),&vals);CHKERRQ(ierr);
118590ace30eSBarry Smith 
118690ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
11870752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
118890ace30eSBarry Smith 
118990ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
119090ace30eSBarry Smith     vals_ptr = vals;
119190ace30eSBarry Smith     for (i=0; i<m; i++) {
119290ace30eSBarry Smith       for (j=0; j<N; j++) {
119390ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
119490ace30eSBarry Smith       }
119590ace30eSBarry Smith     }
119690ace30eSBarry Smith 
119790ace30eSBarry Smith     /* read in other processors and ship out */
119890ace30eSBarry Smith     for (i=1; i<size; i++) {
119990ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12000752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1201ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
120290ace30eSBarry Smith     }
12033a40ed3dSBarry Smith   } else {
120490ace30eSBarry Smith     /* receive numeric values */
1205b0a32e0cSBarry Smith     ierr = PetscMalloc(m*N*sizeof(Scalar),&vals);CHKERRQ(ierr);
120690ace30eSBarry Smith 
120790ace30eSBarry Smith     /* receive message of values*/
1208ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
120990ace30eSBarry Smith 
121090ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
121190ace30eSBarry Smith     vals_ptr = vals;
121290ace30eSBarry Smith     for (i=0; i<m; i++) {
121390ace30eSBarry Smith       for (j=0; j<N; j++) {
121490ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
121590ace30eSBarry Smith       }
121690ace30eSBarry Smith     }
121790ace30eSBarry Smith   }
1218606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1219606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
12206d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12216d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12223a40ed3dSBarry Smith   PetscFunctionReturn(0);
122390ace30eSBarry Smith }
122490ace30eSBarry Smith 
1225273d9f13SBarry Smith EXTERN_C_BEGIN
12265615d1e5SSatish Balay #undef __FUNC__
1227b0a32e0cSBarry Smith #define __FUNC__ "MatLoad_MPIDense"
1228b0a32e0cSBarry Smith int MatLoad_MPIDense(PetscViewer viewer,MatType type,Mat *newmat)
12298965ea79SLois Curfman McInnes {
12308965ea79SLois Curfman McInnes   Mat          A;
12318965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
123219bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
12338965ea79SLois Curfman McInnes   MPI_Status   status;
12348965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
12358965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols;
123619bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
12373a40ed3dSBarry Smith   int          i,nz,ierr,j,rstart,rend,fd;
12388965ea79SLois Curfman McInnes 
12393a40ed3dSBarry Smith   PetscFunctionBegin;
1240d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1241d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12428965ea79SLois Curfman McInnes   if (!rank) {
1243b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
12440752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
124529bbc08cSBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
12468965ea79SLois Curfman McInnes   }
12478965ea79SLois Curfman McInnes 
1248ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
124990ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
125090ace30eSBarry Smith 
125190ace30eSBarry Smith   /*
125290ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
125390ace30eSBarry Smith   */
125490ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
12553a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
12563a40ed3dSBarry Smith     PetscFunctionReturn(0);
125790ace30eSBarry Smith   }
125890ace30eSBarry Smith 
12598965ea79SLois Curfman McInnes   /* determine ownership of all rows */
12608965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
1261b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1262ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
12638965ea79SLois Curfman McInnes   rowners[0] = 0;
12648965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
12658965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
12668965ea79SLois Curfman McInnes   }
12678965ea79SLois Curfman McInnes   rstart = rowners[rank];
12688965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
12698965ea79SLois Curfman McInnes 
12708965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
1271b0a32e0cSBarry Smith   ierr    = PetscMalloc(2*(rend-rstart)*sizeof(int),&ourlens);CHKERRQ(ierr);
12728965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
12738965ea79SLois Curfman McInnes   if (!rank) {
1274b0a32e0cSBarry Smith     ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr);
12750752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1276b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr);
12778965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1278ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1279606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1280ca161407SBarry Smith   } else {
1281ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
12828965ea79SLois Curfman McInnes   }
12838965ea79SLois Curfman McInnes 
12848965ea79SLois Curfman McInnes   if (!rank) {
12858965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
1286b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr);
1287549d3d68SSatish Balay     ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
12888965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
12898965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
12908965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
12918965ea79SLois Curfman McInnes       }
12928965ea79SLois Curfman McInnes     }
1293606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
12948965ea79SLois Curfman McInnes 
12958965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
12968965ea79SLois Curfman McInnes     maxnz = 0;
12978965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
12980452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
12998965ea79SLois Curfman McInnes     }
1300b0a32e0cSBarry Smith     ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr);
13018965ea79SLois Curfman McInnes 
13028965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13038965ea79SLois Curfman McInnes     nz = procsnz[0];
1304b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13050752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13068965ea79SLois Curfman McInnes 
13078965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13088965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13098965ea79SLois Curfman McInnes       nz   = procsnz[i];
13100752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1311ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13128965ea79SLois Curfman McInnes     }
1313606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13143a40ed3dSBarry Smith   } else {
13158965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13168965ea79SLois Curfman McInnes     nz = 0;
13178965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13188965ea79SLois Curfman McInnes       nz += ourlens[i];
13198965ea79SLois Curfman McInnes     }
1320b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13218965ea79SLois Curfman McInnes 
13228965ea79SLois Curfman McInnes     /* receive message of column indices*/
1323ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1324ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
132529bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13268965ea79SLois Curfman McInnes   }
13278965ea79SLois Curfman McInnes 
13288965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1329549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
13308965ea79SLois Curfman McInnes   jj = 0;
13318965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13328965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
13338965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
13348965ea79SLois Curfman McInnes       jj++;
13358965ea79SLois Curfman McInnes     }
13368965ea79SLois Curfman McInnes   }
13378965ea79SLois Curfman McInnes 
13388965ea79SLois Curfman McInnes   /* create our matrix */
13398965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13408965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
13418965ea79SLois Curfman McInnes   }
1342b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
13438965ea79SLois Curfman McInnes   A = *newmat;
13448965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13458965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
13468965ea79SLois Curfman McInnes   }
13478965ea79SLois Curfman McInnes 
13488965ea79SLois Curfman McInnes   if (!rank) {
1349b0a32e0cSBarry Smith     ierr = PetscMalloc(maxnz*sizeof(Scalar),&vals);CHKERRQ(ierr);
13508965ea79SLois Curfman McInnes 
13518965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
13528965ea79SLois Curfman McInnes     nz = procsnz[0];
13530752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
13548965ea79SLois Curfman McInnes 
13558965ea79SLois Curfman McInnes     /* insert into matrix */
13568965ea79SLois Curfman McInnes     jj      = rstart;
13578965ea79SLois Curfman McInnes     smycols = mycols;
13588965ea79SLois Curfman McInnes     svals   = vals;
13598965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13608965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13618965ea79SLois Curfman McInnes       smycols += ourlens[i];
13628965ea79SLois Curfman McInnes       svals   += ourlens[i];
13638965ea79SLois Curfman McInnes       jj++;
13648965ea79SLois Curfman McInnes     }
13658965ea79SLois Curfman McInnes 
13668965ea79SLois Curfman McInnes     /* read in other processors and ship out */
13678965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13688965ea79SLois Curfman McInnes       nz   = procsnz[i];
13690752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1370ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
13718965ea79SLois Curfman McInnes     }
1372606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
13733a40ed3dSBarry Smith   } else {
13748965ea79SLois Curfman McInnes     /* receive numeric values */
1375b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(Scalar),&vals);CHKERRQ(ierr);
13768965ea79SLois Curfman McInnes 
13778965ea79SLois Curfman McInnes     /* receive message of values*/
1378ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1379ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
138029bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13818965ea79SLois Curfman McInnes 
13828965ea79SLois Curfman McInnes     /* insert into matrix */
13838965ea79SLois Curfman McInnes     jj      = rstart;
13848965ea79SLois Curfman McInnes     smycols = mycols;
13858965ea79SLois Curfman McInnes     svals   = vals;
13868965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13878965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13888965ea79SLois Curfman McInnes       smycols += ourlens[i];
13898965ea79SLois Curfman McInnes       svals   += ourlens[i];
13908965ea79SLois Curfman McInnes       jj++;
13918965ea79SLois Curfman McInnes     }
13928965ea79SLois Curfman McInnes   }
1393606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1394606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1395606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1396606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
13978965ea79SLois Curfman McInnes 
13986d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13996d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14003a40ed3dSBarry Smith   PetscFunctionReturn(0);
14018965ea79SLois Curfman McInnes }
1402273d9f13SBarry Smith EXTERN_C_END
140390ace30eSBarry Smith 
140490ace30eSBarry Smith 
140590ace30eSBarry Smith 
140690ace30eSBarry Smith 
140790ace30eSBarry Smith 
1408