xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision c38d4ed214221df9ea04de46f7761bef149d00ff)
1*c38d4ed2SBarry Smith /*$Id: mpidense.c,v 1.131 1999/11/05 14:45:15 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);
199*c38d4ed2SBarry Smith   PLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs);
2003a40ed3dSBarry Smith   PetscFunctionReturn(0);
2018965ea79SLois Curfman McInnes }
2028965ea79SLois Curfman McInnes 
2035615d1e5SSatish Balay #undef __FUNC__
2045615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIDense"
2058f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2068965ea79SLois Curfman McInnes {
20739ddd567SLois Curfman McInnes   Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
2087ef1d9bdSSatish Balay   int          i,n,ierr,*row,*col,flg,j,rstart,ncols;
2097ef1d9bdSSatish Balay   Scalar       *val;
210e0fa3b82SLois Curfman McInnes   InsertMode   addv=mat->insertmode;
2118965ea79SLois Curfman McInnes 
2123a40ed3dSBarry Smith   PetscFunctionBegin;
2138965ea79SLois Curfman McInnes   /*  wait on receives */
2147ef1d9bdSSatish Balay   while (1) {
2158798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2167ef1d9bdSSatish Balay     if (!flg) break;
2178965ea79SLois Curfman McInnes 
2187ef1d9bdSSatish Balay     for ( i=0; i<n; ) {
2197ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2207ef1d9bdSSatish Balay       for ( j=i,rstart=row[j]; j<n; j++ ) { if (row[j] != rstart) break; }
2217ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2227ef1d9bdSSatish Balay       else       ncols = n-i;
2237ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2247ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2257ef1d9bdSSatish Balay       i = j;
2268965ea79SLois Curfman McInnes     }
2277ef1d9bdSSatish Balay   }
2288798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2298965ea79SLois Curfman McInnes 
23039ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
23139ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2328965ea79SLois Curfman McInnes 
2336d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
23439ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2358965ea79SLois Curfman McInnes   }
2363a40ed3dSBarry Smith   PetscFunctionReturn(0);
2378965ea79SLois Curfman McInnes }
2388965ea79SLois Curfman McInnes 
2395615d1e5SSatish Balay #undef __FUNC__
2405615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIDense"
2418f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A)
2428965ea79SLois Curfman McInnes {
2433a40ed3dSBarry Smith   int          ierr;
24439ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense *) A->data;
2453a40ed3dSBarry Smith 
2463a40ed3dSBarry Smith   PetscFunctionBegin;
2473a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2483a40ed3dSBarry Smith   PetscFunctionReturn(0);
2498965ea79SLois Curfman McInnes }
2508965ea79SLois Curfman McInnes 
2515615d1e5SSatish Balay #undef __FUNC__
2525615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_MPIDense"
2538f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs)
2544e220ebcSLois Curfman McInnes {
2553a40ed3dSBarry Smith   PetscFunctionBegin;
2564e220ebcSLois Curfman McInnes   *bs = 1;
2573a40ed3dSBarry Smith   PetscFunctionReturn(0);
2584e220ebcSLois Curfman McInnes }
2594e220ebcSLois Curfman McInnes 
2608965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2618965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2628965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2633501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2648965ea79SLois Curfman McInnes    routine.
2658965ea79SLois Curfman McInnes */
2665615d1e5SSatish Balay #undef __FUNC__
2675615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIDense"
2688f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag)
2698965ea79SLois Curfman McInnes {
27039ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense *) A->data;
2718965ea79SLois Curfman McInnes   int            i,ierr,N, *rows,*owners = l->rowners,size = l->size;
2728965ea79SLois Curfman McInnes   int            *procs,*nprocs,j,found,idx,nsends,*work;
2738965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2748965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2758965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2768965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2778965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2788965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2798965ea79SLois Curfman McInnes   IS             istmp;
2808965ea79SLois Curfman McInnes 
2813a40ed3dSBarry Smith   PetscFunctionBegin;
28277c4ece6SBarry Smith   ierr = ISGetSize(is,&N);CHKERRQ(ierr);
2838965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
2848965ea79SLois Curfman McInnes 
2858965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
2860452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) );CHKPTRQ(nprocs);
287549d3d68SSatish Balay   ierr   = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
288549d3d68SSatish Balay   procs  = nprocs + size;
2890452661fSBarry Smith   owner  = (int *) PetscMalloc((N+1)*sizeof(int));CHKPTRQ(owner); /* see note*/
2908965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
2918965ea79SLois Curfman McInnes     idx = rows[i];
2928965ea79SLois Curfman McInnes     found = 0;
2938965ea79SLois Curfman McInnes     for ( j=0; j<size; j++ ) {
2948965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
2958965ea79SLois Curfman McInnes         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
2968965ea79SLois Curfman McInnes       }
2978965ea79SLois Curfman McInnes     }
298a8c6a408SBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range");
2998965ea79SLois Curfman McInnes   }
3008965ea79SLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
3018965ea79SLois Curfman McInnes 
3028965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
3036831982aSBarry Smith   work   = (int *) PetscMalloc( 2*size*sizeof(int) );CHKPTRQ(work);
3046831982aSBarry Smith   ierr   = MPI_Allreduce( nprocs, work,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr);
3058965ea79SLois Curfman McInnes   nmax   = work[rank];
3066831982aSBarry Smith   nrecvs = work[size+rank];
307606d414cSSatish Balay   ierr   = PetscFree(work);CHKERRQ(ierr);
3088965ea79SLois Curfman McInnes 
3098965ea79SLois Curfman McInnes   /* post receives:   */
3103a40ed3dSBarry Smith   rvalues    = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues);
3113a40ed3dSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits);
3128965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
313ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3148965ea79SLois Curfman McInnes   }
3158965ea79SLois Curfman McInnes 
3168965ea79SLois Curfman McInnes   /* do sends:
3178965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3188965ea79SLois Curfman McInnes          the ith processor
3198965ea79SLois Curfman McInnes   */
3200452661fSBarry Smith   svalues    = (int *) PetscMalloc( (N+1)*sizeof(int) );CHKPTRQ(svalues);
3217056b6fcSBarry Smith   send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
3220452661fSBarry Smith   starts     = (int *) PetscMalloc( (size+1)*sizeof(int) );CHKPTRQ(starts);
3238965ea79SLois Curfman McInnes   starts[0]  = 0;
3248965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3258965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
3268965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3278965ea79SLois Curfman McInnes   }
3288965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3298965ea79SLois Curfman McInnes 
3308965ea79SLois Curfman McInnes   starts[0] = 0;
3318965ea79SLois Curfman McInnes   for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3328965ea79SLois Curfman McInnes   count = 0;
3338965ea79SLois Curfman McInnes   for ( i=0; i<size; i++ ) {
3348965ea79SLois Curfman McInnes     if (procs[i]) {
335ca161407SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3368965ea79SLois Curfman McInnes     }
3378965ea79SLois Curfman McInnes   }
338606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3398965ea79SLois Curfman McInnes 
3408965ea79SLois Curfman McInnes   base = owners[rank];
3418965ea79SLois Curfman McInnes 
3428965ea79SLois Curfman McInnes   /*  wait on receives */
3430452661fSBarry Smith   lens   = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) );CHKPTRQ(lens);
3448965ea79SLois Curfman McInnes   source = lens + nrecvs;
3458965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3468965ea79SLois Curfman McInnes   while (count) {
347ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3488965ea79SLois Curfman McInnes     /* unpack receives into our local space */
349ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3508965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3518965ea79SLois Curfman McInnes     lens[imdex]  = n;
3528965ea79SLois Curfman McInnes     slen += n;
3538965ea79SLois Curfman McInnes     count--;
3548965ea79SLois Curfman McInnes   }
355606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3568965ea79SLois Curfman McInnes 
3578965ea79SLois Curfman McInnes   /* move the data into the send scatter */
3580452661fSBarry Smith   lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) );CHKPTRQ(lrows);
3598965ea79SLois Curfman McInnes   count = 0;
3608965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
3618965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3628965ea79SLois Curfman McInnes     for ( j=0; j<lens[i]; j++ ) {
3638965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3648965ea79SLois Curfman McInnes     }
3658965ea79SLois Curfman McInnes   }
366606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
367606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
368606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
369606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3708965ea79SLois Curfman McInnes 
3718965ea79SLois Curfman McInnes   /* actually zap the local rows */
372029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
3738965ea79SLois Curfman McInnes   PLogObjectParent(A,istmp);
374606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3758965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3768965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes 
3788965ea79SLois Curfman McInnes   /* wait on sends */
3798965ea79SLois Curfman McInnes   if (nsends) {
3807056b6fcSBarry Smith     send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
381ca161407SBarry Smith     ierr        = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
382606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
3838965ea79SLois Curfman McInnes   }
384606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
385606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
3868965ea79SLois Curfman McInnes 
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
3888965ea79SLois Curfman McInnes }
3898965ea79SLois Curfman McInnes 
3905615d1e5SSatish Balay #undef __FUNC__
3915615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIDense"
3928f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3938965ea79SLois Curfman McInnes {
39439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
3958965ea79SLois Curfman McInnes   int          ierr;
396c456f294SBarry Smith 
3973a40ed3dSBarry Smith   PetscFunctionBegin;
39843a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39943a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40044cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4013a40ed3dSBarry Smith   PetscFunctionReturn(0);
4028965ea79SLois Curfman McInnes }
4038965ea79SLois Curfman McInnes 
4045615d1e5SSatish Balay #undef __FUNC__
4055615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIDense"
4068f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4078965ea79SLois Curfman McInnes {
40839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4098965ea79SLois Curfman McInnes   int          ierr;
410c456f294SBarry Smith 
4113a40ed3dSBarry Smith   PetscFunctionBegin;
41243a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41343a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41444cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4153a40ed3dSBarry Smith   PetscFunctionReturn(0);
4168965ea79SLois Curfman McInnes }
4178965ea79SLois Curfman McInnes 
4185615d1e5SSatish Balay #undef __FUNC__
4195615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIDense"
4208f6be9afSLois Curfman McInnes int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy)
421096963f5SLois Curfman McInnes {
422096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
423096963f5SLois Curfman McInnes   int          ierr;
4243501a2bdSLois Curfman McInnes   Scalar       zero = 0.0;
425096963f5SLois Curfman McInnes 
4263a40ed3dSBarry Smith   PetscFunctionBegin;
4273501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
42844cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
429537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
430537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4313a40ed3dSBarry Smith   PetscFunctionReturn(0);
432096963f5SLois Curfman McInnes }
433096963f5SLois Curfman McInnes 
4345615d1e5SSatish Balay #undef __FUNC__
4355615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIDense"
4368f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
437096963f5SLois Curfman McInnes {
438096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
439096963f5SLois Curfman McInnes   int          ierr;
440096963f5SLois Curfman McInnes 
4413a40ed3dSBarry Smith   PetscFunctionBegin;
4423501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
44344cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
444537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
445537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4463a40ed3dSBarry Smith   PetscFunctionReturn(0);
447096963f5SLois Curfman McInnes }
448096963f5SLois Curfman McInnes 
4495615d1e5SSatish Balay #undef __FUNC__
4505615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIDense"
4518f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v)
4528965ea79SLois Curfman McInnes {
45339ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
454096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data;
45544cd7ae7SLois Curfman McInnes   int          ierr, len, i, n, m = a->m, radd;
45644cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
457ed3cc1f0SBarry Smith 
4583a40ed3dSBarry Smith   PetscFunctionBegin;
45944cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
460096963f5SLois Curfman McInnes   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
461096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
462a8c6a408SBarry Smith   if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
46344cd7ae7SLois Curfman McInnes   len = PetscMin(aloc->m,aloc->n);
4647ddc982cSLois Curfman McInnes   radd = a->rstart*m;
46544cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
466096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
467096963f5SLois Curfman McInnes   }
4689a8c540fSSatish Balay   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
4693a40ed3dSBarry Smith   PetscFunctionReturn(0);
4708965ea79SLois Curfman McInnes }
4718965ea79SLois Curfman McInnes 
4725615d1e5SSatish Balay #undef __FUNC__
4735615d1e5SSatish Balay #define __FUNC__ "MatDestroy_MPIDense"
474e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat)
4758965ea79SLois Curfman McInnes {
4763501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4778965ea79SLois Curfman McInnes   int          ierr;
478ed3cc1f0SBarry Smith 
4793a40ed3dSBarry Smith   PetscFunctionBegin;
48094d884c6SBarry Smith 
48194d884c6SBarry Smith   if (mat->mapping) {
48294d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr);
48394d884c6SBarry Smith   }
48494d884c6SBarry Smith   if (mat->bmapping) {
48594d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr);
48694d884c6SBarry Smith   }
487aa482453SBarry Smith #if defined(PETSC_USE_LOG)
488e1311b90SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N);
4898965ea79SLois Curfman McInnes #endif
4908798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
491606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
4923501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
4933501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4943501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
495622d7880SLois Curfman McInnes   if (mdn->factor) {
496606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
497606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
498606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
499606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
500622d7880SLois Curfman McInnes   }
501606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
50261b13de0SBarry Smith   if (mat->rmap) {
50361b13de0SBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
50461b13de0SBarry Smith   }
50561b13de0SBarry Smith   if (mat->cmap) {
50661b13de0SBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
50790f02eecSBarry Smith   }
5088965ea79SLois Curfman McInnes   PLogObjectDestroy(mat);
5090452661fSBarry Smith   PetscHeaderDestroy(mat);
5103a40ed3dSBarry Smith   PetscFunctionReturn(0);
5118965ea79SLois Curfman McInnes }
51239ddd567SLois Curfman McInnes 
5135615d1e5SSatish Balay #undef __FUNC__
5145615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary"
51539ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer)
5168965ea79SLois Curfman McInnes {
51739ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
5188965ea79SLois Curfman McInnes   int          ierr;
5197056b6fcSBarry Smith 
5203a40ed3dSBarry Smith   PetscFunctionBegin;
52139ddd567SLois Curfman McInnes   if (mdn->size == 1) {
52239ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5238965ea79SLois Curfman McInnes   }
524a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported");
5253a40ed3dSBarry Smith   PetscFunctionReturn(0);
5268965ea79SLois Curfman McInnes }
5278965ea79SLois Curfman McInnes 
5285615d1e5SSatish Balay #undef __FUNC__
529f1af5d2fSBarry Smith #define __FUNC__ "MatView_MPIDense_ASCIIorDraworSocket"
530f1af5d2fSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,Viewer viewer)
5318965ea79SLois Curfman McInnes {
53239ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
53377ed5343SBarry Smith   int          ierr, format, size = mdn->size, rank = mdn->rank;
53419bcc07fSBarry Smith   ViewerType   vtype;
535f1af5d2fSBarry Smith   PetscTruth   isascii,isdraw;
5368965ea79SLois Curfman McInnes 
5373a40ed3dSBarry Smith   PetscFunctionBegin;
538f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
539f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
540f1af5d2fSBarry Smith   if (isascii) {
5413a40ed3dSBarry Smith     ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
542888f2ed8SSatish Balay     ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
543639f9d9dSBarry Smith     if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
5444e220ebcSLois Curfman McInnes       MatInfo info;
545888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
5466831982aSBarry Smith       ierr = ViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m,
5476831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
5486831982aSBarry Smith       ierr = ViewerFlush(viewer);CHKERRQ(ierr);
5493501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5503a40ed3dSBarry Smith       PetscFunctionReturn(0);
55196f6c058SBarry Smith     } else if (format == VIEWER_FORMAT_ASCII_INFO) {
5523a40ed3dSBarry Smith       PetscFunctionReturn(0);
5538965ea79SLois Curfman McInnes     }
554f1af5d2fSBarry Smith   } else if (isdraw) {
555f1af5d2fSBarry Smith     Draw       draw;
556f1af5d2fSBarry Smith     PetscTruth isnull;
557f1af5d2fSBarry Smith 
558f1af5d2fSBarry Smith     ierr = ViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
559f1af5d2fSBarry Smith     ierr = DrawIsNull(draw, &isnull);CHKERRQ(ierr);
560f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
561f1af5d2fSBarry Smith   }
56277ed5343SBarry Smith 
5638965ea79SLois Curfman McInnes   if (size == 1) {
56439ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5653a40ed3dSBarry Smith   } else {
5668965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5678965ea79SLois Curfman McInnes     Mat          A;
56839ddd567SLois Curfman McInnes     int          M = mdn->M, N = mdn->N,m,row,i, nz, *cols;
56939ddd567SLois Curfman McInnes     Scalar       *vals;
57039ddd567SLois Curfman McInnes     Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data;
5718965ea79SLois Curfman McInnes 
5728965ea79SLois Curfman McInnes     if (!rank) {
573f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5743a40ed3dSBarry Smith     } else {
575f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5768965ea79SLois Curfman McInnes     }
5778965ea79SLois Curfman McInnes     PLogObjectParent(mat,A);
5788965ea79SLois Curfman McInnes 
57939ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
58039ddd567SLois Curfman McInnes        but it's quick for now */
58139ddd567SLois Curfman McInnes     row = mdn->rstart; m = Amdn->m;
5828965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
58339ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
58439ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
58539ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
58639ddd567SLois Curfman McInnes       row++;
5878965ea79SLois Curfman McInnes     }
5888965ea79SLois Curfman McInnes 
5896d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5906d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5918965ea79SLois Curfman McInnes     if (!rank) {
5926831982aSBarry Smith       Viewer sviewer;
5936831982aSBarry Smith       ierr = ViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
5946831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5956831982aSBarry Smith       ierr = ViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
5968965ea79SLois Curfman McInnes     }
5976831982aSBarry Smith     ierr = ViewerFlush(viewer);CHKERRQ(ierr);
5988965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
5998965ea79SLois Curfman McInnes   }
6003a40ed3dSBarry Smith   PetscFunctionReturn(0);
6018965ea79SLois Curfman McInnes }
6028965ea79SLois Curfman McInnes 
6035615d1e5SSatish Balay #undef __FUNC__
6045615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense"
605e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer)
6068965ea79SLois Curfman McInnes {
60739ddd567SLois Curfman McInnes   int        ierr;
608f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
6098965ea79SLois Curfman McInnes 
610433994e6SBarry Smith   PetscFunctionBegin;
6110f5bd95cSBarry Smith 
6126831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
6136831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
614f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
615f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
6160f5bd95cSBarry Smith 
617f1af5d2fSBarry Smith   if (isascii || issocket || isdraw) {
618f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6190f5bd95cSBarry Smith   } else if (isbinary) {
6203a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6215cd90555SBarry Smith   } else {
6220f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6238965ea79SLois Curfman McInnes   }
6243a40ed3dSBarry Smith   PetscFunctionReturn(0);
6258965ea79SLois Curfman McInnes }
6268965ea79SLois Curfman McInnes 
6275615d1e5SSatish Balay #undef __FUNC__
6285615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense"
6298f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6308965ea79SLois Curfman McInnes {
6313501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6323501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6334e220ebcSLois Curfman McInnes   int          ierr;
6344e220ebcSLois Curfman McInnes   double       isend[5], irecv[5];
6358965ea79SLois Curfman McInnes 
6363a40ed3dSBarry Smith   PetscFunctionBegin;
6374e220ebcSLois Curfman McInnes   info->rows_global    = (double)mat->M;
6384e220ebcSLois Curfman McInnes   info->columns_global = (double)mat->N;
6394e220ebcSLois Curfman McInnes   info->rows_local     = (double)mat->m;
6404e220ebcSLois Curfman McInnes   info->columns_local  = (double)mat->N;
6414e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6424e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6434e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6444e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6458965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6464e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6474e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6484e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6494e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6504e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6518965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
652f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
6534e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6544e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6554e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6564e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6574e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6588965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
659f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
6604e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6614e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6624e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6634e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6644e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6658965ea79SLois Curfman McInnes   }
6664e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6674e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6684e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6693a40ed3dSBarry Smith   PetscFunctionReturn(0);
6708965ea79SLois Curfman McInnes }
6718965ea79SLois Curfman McInnes 
6728c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*);
6738aaee692SLois Curfman McInnes    extern int MatLUFactorNumeric_MPIDense(Mat,Mat*);
6748aaee692SLois Curfman McInnes    extern int MatLUFactor_MPIDense(Mat,IS,IS,double);
6758aaee692SLois Curfman McInnes    extern int MatSolve_MPIDense(Mat,Vec,Vec);
6768c469469SLois Curfman McInnes    extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec);
6778aaee692SLois Curfman McInnes    extern int MatSolveTrans_MPIDense(Mat,Vec,Vec);
6788aaee692SLois Curfman McInnes    extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */
6798aaee692SLois Curfman McInnes 
6805615d1e5SSatish Balay #undef __FUNC__
6815615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense"
6828f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6838965ea79SLois Curfman McInnes {
68439ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
6858965ea79SLois Curfman McInnes 
6863a40ed3dSBarry Smith   PetscFunctionBegin;
6876d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6886d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
6894787f768SSatish Balay       op == MAT_NEW_NONZERO_LOCATION_ERR ||
6904787f768SSatish Balay       op == MAT_NEW_NONZERO_ALLOCATION_ERR ||
691219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
692219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED) {
693b1fbbac0SLois Curfman McInnes         MatSetOption(a->A,op);
694b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
695aeafbbfcSLois Curfman McInnes         a->roworiented = 1;
6968965ea79SLois Curfman McInnes         MatSetOption(a->A,op);
697b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
698219d9a1aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
6996d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
7006d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
701b51ba29fSSatish Balay              op == MAT_YES_NEW_DIAGONALS ||
702b51ba29fSSatish Balay              op == MAT_USE_HASH_TABLE) {
703981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
7043a40ed3dSBarry Smith   } else if (op == MAT_COLUMN_ORIENTED) {
7053a40ed3dSBarry Smith     a->roworiented = 0; MatSetOption(a->A,op);
7063782ba37SSatish Balay   } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) {
7073782ba37SSatish Balay     a->donotstash = 1;
7083a40ed3dSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS) {
7093a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
7103a40ed3dSBarry Smith   } else {
7113a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
7123a40ed3dSBarry Smith   }
7133a40ed3dSBarry Smith   PetscFunctionReturn(0);
7148965ea79SLois Curfman McInnes }
7158965ea79SLois Curfman McInnes 
7165615d1e5SSatish Balay #undef __FUNC__
7175615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense"
7188f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n)
7198965ea79SLois Curfman McInnes {
7203501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7213a40ed3dSBarry Smith 
7223a40ed3dSBarry Smith   PetscFunctionBegin;
723f1af5d2fSBarry Smith   if (m) *m = mat->M;
724f1af5d2fSBarry Smith   if (n) *n = mat->N;
7253a40ed3dSBarry Smith   PetscFunctionReturn(0);
7268965ea79SLois Curfman McInnes }
7278965ea79SLois Curfman McInnes 
7285615d1e5SSatish Balay #undef __FUNC__
7295615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense"
7308f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n)
7318965ea79SLois Curfman McInnes {
7323501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7333a40ed3dSBarry Smith 
7343a40ed3dSBarry Smith   PetscFunctionBegin;
7358965ea79SLois Curfman McInnes   *m = mat->m; *n = mat->N;
7363a40ed3dSBarry Smith   PetscFunctionReturn(0);
7378965ea79SLois Curfman McInnes }
7388965ea79SLois Curfman McInnes 
7395615d1e5SSatish Balay #undef __FUNC__
7405615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense"
7418f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
7428965ea79SLois Curfman McInnes {
7433501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7443a40ed3dSBarry Smith 
7453a40ed3dSBarry Smith   PetscFunctionBegin;
7468965ea79SLois Curfman McInnes   *m = mat->rstart; *n = mat->rend;
7473a40ed3dSBarry Smith   PetscFunctionReturn(0);
7488965ea79SLois Curfman McInnes }
7498965ea79SLois Curfman McInnes 
7505615d1e5SSatish Balay #undef __FUNC__
7515615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense"
7528f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
7538965ea79SLois Curfman McInnes {
7543501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7553a40ed3dSBarry Smith   int          lrow, rstart = mat->rstart, rend = mat->rend,ierr;
7568965ea79SLois Curfman McInnes 
7573a40ed3dSBarry Smith   PetscFunctionBegin;
758a8c6a408SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows")
7598965ea79SLois Curfman McInnes   lrow = row - rstart;
7603a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7613a40ed3dSBarry Smith   PetscFunctionReturn(0);
7628965ea79SLois Curfman McInnes }
7638965ea79SLois Curfman McInnes 
7645615d1e5SSatish Balay #undef __FUNC__
7655615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense"
7668f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
7678965ea79SLois Curfman McInnes {
768606d414cSSatish Balay   int ierr;
769606d414cSSatish Balay 
7703a40ed3dSBarry Smith   PetscFunctionBegin;
771606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
772606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7733a40ed3dSBarry Smith   PetscFunctionReturn(0);
7748965ea79SLois Curfman McInnes }
7758965ea79SLois Curfman McInnes 
7765615d1e5SSatish Balay #undef __FUNC__
7775b2fa520SLois Curfman McInnes #define __FUNC__ "MatDiagonalScale_MPIDense"
7785b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7795b2fa520SLois Curfman McInnes {
7805b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
7815b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data;
7825b2fa520SLois Curfman McInnes   Scalar       *l,*r,x,*v;
78372d926a5SLois Curfman McInnes   int          ierr,i,j,s2a,s3a,s2,s3,m=mat->m,n=mat->n;
7845b2fa520SLois Curfman McInnes 
7855b2fa520SLois Curfman McInnes   PetscFunctionBegin;
78672d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7875b2fa520SLois Curfman McInnes   if (ll) {
78872d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
78972d926a5SLois Curfman McInnes     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector non-conforming local size");
7905b2fa520SLois Curfman McInnes     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7915b2fa520SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
7925b2fa520SLois Curfman McInnes       x = l[i];
7935b2fa520SLois Curfman McInnes       v = mat->v + i;
7945b2fa520SLois Curfman McInnes       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
7955b2fa520SLois Curfman McInnes     }
7965b2fa520SLois Curfman McInnes     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
7975b2fa520SLois Curfman McInnes     PLogFlops(n*m);
7985b2fa520SLois Curfman McInnes   }
7995b2fa520SLois Curfman McInnes   if (rr) {
80072d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
80172d926a5SLois Curfman McInnes     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec non-conforming local size");
8025b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8035b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8045b2fa520SLois Curfman McInnes     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8055b2fa520SLois Curfman McInnes     for ( i=0; i<n; i++ ) {
8065b2fa520SLois Curfman McInnes       x = r[i];
8075b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8085b2fa520SLois Curfman McInnes       for ( j=0; j<m; j++ ) { (*v++) *= x;}
8095b2fa520SLois Curfman McInnes     }
81072d926a5SLois Curfman McInnes     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
8115b2fa520SLois Curfman McInnes     PLogFlops(n*m);
8125b2fa520SLois Curfman McInnes   }
8135b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8145b2fa520SLois Curfman McInnes }
8155b2fa520SLois Curfman McInnes 
8165b2fa520SLois Curfman McInnes #undef __FUNC__
8175615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense"
8188f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm)
819096963f5SLois Curfman McInnes {
8203501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
8213501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data;
8223501a2bdSLois Curfman McInnes   int          ierr, i, j;
8233501a2bdSLois Curfman McInnes   double       sum = 0.0;
8243501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
8253501a2bdSLois Curfman McInnes 
8263a40ed3dSBarry Smith   PetscFunctionBegin;
8273501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
8283501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm);CHKERRQ(ierr);
8293501a2bdSLois Curfman McInnes   } else {
8303501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
8313501a2bdSLois Curfman McInnes       for (i=0; i<mat->n*mat->m; i++ ) {
832aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
833e20fef11SSatish Balay         sum += PetscReal(PetscConj(*v)*(*v)); v++;
8343501a2bdSLois Curfman McInnes #else
8353501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
8363501a2bdSLois Curfman McInnes #endif
8373501a2bdSLois Curfman McInnes       }
838ca161407SBarry Smith       ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
8393501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
8403501a2bdSLois Curfman McInnes       PLogFlops(2*mat->n*mat->m);
8413a40ed3dSBarry Smith     } else if (type == NORM_1) {
8423501a2bdSLois Curfman McInnes       double *tmp, *tmp2;
8430452661fSBarry Smith       tmp  = (double *) PetscMalloc( 2*mdn->N*sizeof(double) );CHKPTRQ(tmp);
8443501a2bdSLois Curfman McInnes       tmp2 = tmp + mdn->N;
845549d3d68SSatish Balay       ierr = PetscMemzero(tmp,2*mdn->N*sizeof(double));CHKERRQ(ierr);
846096963f5SLois Curfman McInnes       *norm = 0.0;
8473501a2bdSLois Curfman McInnes       v = mat->v;
8483501a2bdSLois Curfman McInnes       for ( j=0; j<mat->n; j++ ) {
8493501a2bdSLois Curfman McInnes         for ( i=0; i<mat->m; i++ ) {
85067e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8513501a2bdSLois Curfman McInnes         }
8523501a2bdSLois Curfman McInnes       }
853ca161407SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
8543501a2bdSLois Curfman McInnes       for ( j=0; j<mdn->N; j++ ) {
8553501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
8563501a2bdSLois Curfman McInnes       }
857606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
8583501a2bdSLois Curfman McInnes       PLogFlops(mat->n*mat->m);
8593a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
8603501a2bdSLois Curfman McInnes       double ntemp;
8613501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
862ca161407SBarry Smith       ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
8633a40ed3dSBarry Smith     } else {
864a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
8653501a2bdSLois Curfman McInnes     }
8663501a2bdSLois Curfman McInnes   }
8673a40ed3dSBarry Smith   PetscFunctionReturn(0);
8683501a2bdSLois Curfman McInnes }
8693501a2bdSLois Curfman McInnes 
8705615d1e5SSatish Balay #undef __FUNC__
8715615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense"
8728f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8733501a2bdSLois Curfman McInnes {
8743501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
8753501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data;
8763501a2bdSLois Curfman McInnes   Mat          B;
8773501a2bdSLois Curfman McInnes   int          M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart;
8783501a2bdSLois Curfman McInnes   int          j, i, ierr;
8793501a2bdSLois Curfman McInnes   Scalar       *v;
8803501a2bdSLois Curfman McInnes 
8813a40ed3dSBarry Smith   PetscFunctionBegin;
8827056b6fcSBarry Smith   if (matout == PETSC_NULL && M != N) {
883a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place");
8847056b6fcSBarry Smith   }
8857056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8863501a2bdSLois Curfman McInnes 
8873501a2bdSLois Curfman McInnes   m = Aloc->m; n = Aloc->n; v = Aloc->v;
8880452661fSBarry Smith   rwork = (int *) PetscMalloc(n*sizeof(int));CHKPTRQ(rwork);
8893501a2bdSLois Curfman McInnes   for ( j=0; j<n; j++ ) {
8903501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8913501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8923501a2bdSLois Curfman McInnes     v   += m;
8933501a2bdSLois Curfman McInnes   }
894606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8956d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8966d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8973638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
8983501a2bdSLois Curfman McInnes     *matout = B;
8993501a2bdSLois Curfman McInnes   } else {
900f830108cSBarry Smith     PetscOps *Abops;
90109dc0095SBarry Smith     MatOps   Aops;
902f830108cSBarry Smith 
9033501a2bdSLois Curfman McInnes     /* This isn't really an in-place transpose, but free data struct from a */
904606d414cSSatish Balay     ierr = PetscFree(a->rowners);CHKERRQ(ierr);
9053501a2bdSLois Curfman McInnes     ierr = MatDestroy(a->A);CHKERRQ(ierr);
9063501a2bdSLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
9073501a2bdSLois Curfman McInnes     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
908606d414cSSatish Balay     ierr = PetscFree(a);CHKERRQ(ierr);
909f830108cSBarry Smith 
910f830108cSBarry Smith     /*
911f830108cSBarry Smith          This is horrible, horrible code. We need to keep the
912f830108cSBarry Smith       A pointers for the bops and ops but copy everything
913f830108cSBarry Smith       else from C.
914f830108cSBarry Smith     */
915f830108cSBarry Smith     Abops   = A->bops;
916f830108cSBarry Smith     Aops    = A->ops;
917549d3d68SSatish Balay     ierr    = PetscMemcpy(A,B,sizeof(struct _p_Mat));CHKERRQ(ierr);
918f830108cSBarry Smith     A->bops = Abops;
919f830108cSBarry Smith     A->ops  = Aops;
920f830108cSBarry Smith 
9210452661fSBarry Smith     PetscHeaderDestroy(B);
9223501a2bdSLois Curfman McInnes   }
9233a40ed3dSBarry Smith   PetscFunctionReturn(0);
924096963f5SLois Curfman McInnes }
925096963f5SLois Curfman McInnes 
926eadb2fb4SBarry Smith #include "pinclude/blaslapack.h"
9275615d1e5SSatish Balay #undef __FUNC__
9285615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense"
9298f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA)
93044cd7ae7SLois Curfman McInnes {
93144cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense *) inA->data;
93244cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->A->data;
93344cd7ae7SLois Curfman McInnes   int          one = 1, nz;
93444cd7ae7SLois Curfman McInnes 
9353a40ed3dSBarry Smith   PetscFunctionBegin;
93644cd7ae7SLois Curfman McInnes   nz = a->m*a->n;
93744cd7ae7SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
93844cd7ae7SLois Curfman McInnes   PLogFlops(nz);
9393a40ed3dSBarry Smith   PetscFunctionReturn(0);
94044cd7ae7SLois Curfman McInnes }
94144cd7ae7SLois Curfman McInnes 
9425609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
9437b2a1423SBarry Smith extern int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **);
9448965ea79SLois Curfman McInnes 
9458965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
94609dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
94709dc0095SBarry Smith        MatGetRow_MPIDense,
94809dc0095SBarry Smith        MatRestoreRow_MPIDense,
94909dc0095SBarry Smith        MatMult_MPIDense,
95009dc0095SBarry Smith        MatMultAdd_MPIDense,
95109dc0095SBarry Smith        MatMultTrans_MPIDense,
95209dc0095SBarry Smith        MatMultTransAdd_MPIDense,
9538965ea79SLois Curfman McInnes        0,
95409dc0095SBarry Smith        0,
95509dc0095SBarry Smith        0,
95609dc0095SBarry Smith        0,
95709dc0095SBarry Smith        0,
95809dc0095SBarry Smith        0,
95909dc0095SBarry Smith        0,
96009dc0095SBarry Smith        MatTranspose_MPIDense,
96109dc0095SBarry Smith        MatGetInfo_MPIDense,0,
96209dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9635b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
96409dc0095SBarry Smith        MatNorm_MPIDense,
96509dc0095SBarry Smith        MatAssemblyBegin_MPIDense,
96609dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
96709dc0095SBarry Smith        0,
96809dc0095SBarry Smith        MatSetOption_MPIDense,
96909dc0095SBarry Smith        MatZeroEntries_MPIDense,
97009dc0095SBarry Smith        MatZeroRows_MPIDense,
97109dc0095SBarry Smith        0,
97209dc0095SBarry Smith        0,
97309dc0095SBarry Smith        0,
97409dc0095SBarry Smith        0,
97509dc0095SBarry Smith        MatGetSize_MPIDense,
97609dc0095SBarry Smith        MatGetLocalSize_MPIDense,
97739ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
97809dc0095SBarry Smith        0,
97909dc0095SBarry Smith        0,
98009dc0095SBarry Smith        MatGetArray_MPIDense,
98109dc0095SBarry Smith        MatRestoreArray_MPIDense,
9825609ef8eSBarry Smith        MatDuplicate_MPIDense,
98309dc0095SBarry Smith        0,
98409dc0095SBarry Smith        0,
98509dc0095SBarry Smith        0,
98609dc0095SBarry Smith        0,
98709dc0095SBarry Smith        0,
9882ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
98909dc0095SBarry Smith        0,
99009dc0095SBarry Smith        MatGetValues_MPIDense,
99109dc0095SBarry Smith        0,
99209dc0095SBarry Smith        0,
99309dc0095SBarry Smith        MatScale_MPIDense,
99409dc0095SBarry Smith        0,
99509dc0095SBarry Smith        0,
99609dc0095SBarry Smith        0,
99709dc0095SBarry Smith        MatGetBlockSize_MPIDense,
99809dc0095SBarry Smith        0,
99909dc0095SBarry Smith        0,
100009dc0095SBarry Smith        0,
100109dc0095SBarry Smith        0,
100209dc0095SBarry Smith        0,
100309dc0095SBarry Smith        0,
100409dc0095SBarry Smith        0,
100509dc0095SBarry Smith        0,
100609dc0095SBarry Smith        0,
1007ca3fa75bSLois Curfman McInnes        MatGetSubMatrix_MPIDense,
100809dc0095SBarry Smith        0,
100909dc0095SBarry Smith        0,
101009dc0095SBarry Smith        MatGetMaps_Petsc};
10118965ea79SLois Curfman McInnes 
10125615d1e5SSatish Balay #undef __FUNC__
10135615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense"
10148965ea79SLois Curfman McInnes /*@C
101539ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
10168965ea79SLois Curfman McInnes 
1017db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1018db81eaa0SLois Curfman McInnes 
10198965ea79SLois Curfman McInnes    Input Parameters:
1020db81eaa0SLois Curfman McInnes +  comm - MPI communicator
10218965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1022db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
10238965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1024db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1025db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1026dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
10278965ea79SLois Curfman McInnes 
10288965ea79SLois Curfman McInnes    Output Parameter:
1029477f1c0bSLois Curfman McInnes .  A - the matrix
10308965ea79SLois Curfman McInnes 
1031b259b22eSLois Curfman McInnes    Notes:
103239ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
103339ddd567SLois Curfman McInnes    storage by columns.
10348965ea79SLois Curfman McInnes 
103518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
103618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1037b4fd4287SBarry Smith    set data=PETSC_NULL.
103818f449edSLois Curfman McInnes 
10398965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
10408965ea79SLois Curfman McInnes    (possibly both).
10418965ea79SLois Curfman McInnes 
1042027ccd11SLois Curfman McInnes    Level: intermediate
1043027ccd11SLois Curfman McInnes 
104439ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel
10458965ea79SLois Curfman McInnes 
104639ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
10478965ea79SLois Curfman McInnes @*/
1048477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
10498965ea79SLois Curfman McInnes {
10508965ea79SLois Curfman McInnes   Mat          mat;
105139ddd567SLois Curfman McInnes   Mat_MPIDense *a;
1052f1af5d2fSBarry Smith   int          ierr, i;
1053f1af5d2fSBarry Smith   PetscTruth   flg;
10548965ea79SLois Curfman McInnes 
10553a40ed3dSBarry Smith   PetscFunctionBegin;
1056ed2daf61SLois Curfman McInnes   /* Note:  For now, when data is specified above, this assumes the user correctly
1057ed2daf61SLois Curfman McInnes    allocates the local dense storage space.  We should add error checking. */
105818f449edSLois Curfman McInnes 
1059477f1c0bSLois Curfman McInnes   *A = 0;
10603f1db9ecSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",comm,MatDestroy,MatView);
10618965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
10620452661fSBarry Smith   mat->data         = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a);
1063549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1064e1311b90SBarry Smith   mat->ops->destroy = MatDestroy_MPIDense;
1065e1311b90SBarry Smith   mat->ops->view    = MatView_MPIDense;
10668965ea79SLois Curfman McInnes   mat->factor       = 0;
106790f02eecSBarry Smith   mat->mapping      = 0;
10688965ea79SLois Curfman McInnes 
1069622d7880SLois Curfman McInnes   a->factor       = 0;
1070e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1071d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&a->rank);CHKERRQ(ierr);
1072d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&a->size);CHKERRQ(ierr);
10738965ea79SLois Curfman McInnes 
107496f6c058SBarry Smith   ierr = PetscSplitOwnership(comm,&m,&M);CHKERRQ(ierr);
107539ddd567SLois Curfman McInnes 
1076f1af5d2fSBarry Smith   ierr = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr);
1077be0abb6dSBarry Smith   a->nvec = n;
1078c7fcc2eaSBarry Smith 
107939ddd567SLois Curfman McInnes   /* each row stores all columns */
1080aca0ad90SLois Curfman McInnes   a->N = mat->N = N;
1081aca0ad90SLois Curfman McInnes   a->M = mat->M = M;
1082aca0ad90SLois Curfman McInnes   a->m = mat->m = m;
1083be0abb6dSBarry Smith   a->n = mat->n = N;   /* NOTE: n == N */
10848965ea79SLois Curfman McInnes 
1085c7fcc2eaSBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1086c7fcc2eaSBarry Smith      we should remove the duplicate information that is not contained in the maps */
1087488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,M,&mat->rmap);CHKERRQ(ierr);
1088be0abb6dSBarry Smith   ierr = MapCreateMPI(comm,n,N,&mat->cmap);CHKERRQ(ierr);
1089c7fcc2eaSBarry Smith 
10908965ea79SLois Curfman McInnes   /* build local table of row and column ownerships */
1091d7e8b826SBarry Smith   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners);
1092d7e8b826SBarry Smith   a->cowners = a->rowners + a->size + 1;
1093f09e8eb9SSatish Balay   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1094ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
10958965ea79SLois Curfman McInnes   a->rowners[0] = 0;
10968965ea79SLois Curfman McInnes   for ( i=2; i<=a->size; i++ ) {
10978965ea79SLois Curfman McInnes     a->rowners[i] += a->rowners[i-1];
10988965ea79SLois Curfman McInnes   }
10998965ea79SLois Curfman McInnes   a->rstart = a->rowners[a->rank];
11008965ea79SLois Curfman McInnes   a->rend   = a->rowners[a->rank+1];
1101ca161407SBarry Smith   ierr      = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
1102d7e8b826SBarry Smith   a->cowners[0] = 0;
1103d7e8b826SBarry Smith   for ( i=2; i<=a->size; i++ ) {
1104d7e8b826SBarry Smith     a->cowners[i] += a->cowners[i-1];
1105d7e8b826SBarry Smith   }
11068965ea79SLois Curfman McInnes 
1107029af93fSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A);CHKERRQ(ierr);
11088965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
11098965ea79SLois Curfman McInnes 
11108965ea79SLois Curfman McInnes   /* build cache for off array entries formed */
11113782ba37SSatish Balay   a->donotstash = 0;
11128798bf22SSatish Balay   ierr = MatStashCreate_Private(comm,1,&mat->stash);CHKERRQ(ierr);
11138965ea79SLois Curfman McInnes 
11148965ea79SLois Curfman McInnes   /* stuff used for matrix vector multiply */
11158965ea79SLois Curfman McInnes   a->lvec        = 0;
11168965ea79SLois Curfman McInnes   a->Mvctx       = 0;
111739b7565bSBarry Smith   a->roworiented = 1;
11188965ea79SLois Curfman McInnes 
1119f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
11200de54da6SSatish Balay                                      "MatGetDiagonalBlock_MPIDense",
11210de54da6SSatish Balay                                      (void*)MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
11220de54da6SSatish Balay 
1123477f1c0bSLois Curfman McInnes   *A = mat;
112425cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
112525cdf11fSBarry Smith   if (flg) {
11268c469469SLois Curfman McInnes     ierr = MatPrintHelp(mat);CHKERRQ(ierr);
11278c469469SLois Curfman McInnes   }
11283a40ed3dSBarry Smith   PetscFunctionReturn(0);
11298965ea79SLois Curfman McInnes }
11308965ea79SLois Curfman McInnes 
11315615d1e5SSatish Balay #undef __FUNC__
11325609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_MPIDense"
11335609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11348965ea79SLois Curfman McInnes {
11358965ea79SLois Curfman McInnes   Mat          mat;
11363501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data;
113739ddd567SLois Curfman McInnes   int          ierr;
11382ba99913SLois Curfman McInnes   FactorCtx    *factor;
11398965ea79SLois Curfman McInnes 
11403a40ed3dSBarry Smith   PetscFunctionBegin;
11418965ea79SLois Curfman McInnes   *newmat       = 0;
11423f1db9ecSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",A->comm,MatDestroy,MatView);
11438965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
11440452661fSBarry Smith   mat->data         = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a);
1145549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1146e1311b90SBarry Smith   mat->ops->destroy = MatDestroy_MPIDense;
1147e1311b90SBarry Smith   mat->ops->view    = MatView_MPIDense;
11483501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1149c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
11508965ea79SLois Curfman McInnes 
115144cd7ae7SLois Curfman McInnes   a->m = mat->m = oldmat->m;
115244cd7ae7SLois Curfman McInnes   a->n = mat->n = oldmat->n;
115344cd7ae7SLois Curfman McInnes   a->M = mat->M = oldmat->M;
115444cd7ae7SLois Curfman McInnes   a->N = mat->N = oldmat->N;
11552ba99913SLois Curfman McInnes   if (oldmat->factor) {
11562ba99913SLois Curfman McInnes     a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx));CHKPTRQ(factor);
11572ba99913SLois Curfman McInnes     /* copy factor contents ... add this code! */
11582ba99913SLois Curfman McInnes   } else a->factor = 0;
11598965ea79SLois Curfman McInnes 
11608965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
11618965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
11628965ea79SLois Curfman McInnes   a->size         = oldmat->size;
11638965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1164e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
11653782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
11660452661fSBarry Smith   a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners);
1167f09e8eb9SSatish Balay   PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1168549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
11698798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
11708965ea79SLois Curfman McInnes 
11718965ea79SLois Curfman McInnes   ierr =  VecDuplicate(oldmat->lvec,&a->lvec);CHKERRQ(ierr);
11728965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->lvec);
117355659b69SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx);CHKERRQ(ierr);
11748965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->Mvctx);
11755609ef8eSBarry Smith   ierr =  MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
11768965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
11778965ea79SLois Curfman McInnes   *newmat = mat;
11783a40ed3dSBarry Smith   PetscFunctionReturn(0);
11798965ea79SLois Curfman McInnes }
11808965ea79SLois Curfman McInnes 
118177c4ece6SBarry Smith #include "sys.h"
11828965ea79SLois Curfman McInnes 
11835615d1e5SSatish Balay #undef __FUNC__
11845615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile"
118590ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat)
118690ace30eSBarry Smith {
118740011551SBarry Smith   int        *rowners, i,size,rank,m,ierr,nz,j;
118890ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
118990ace30eSBarry Smith   MPI_Status status;
119090ace30eSBarry Smith 
11913a40ed3dSBarry Smith   PetscFunctionBegin;
1192d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1193d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
119490ace30eSBarry Smith 
119590ace30eSBarry Smith   /* determine ownership of all rows */
119690ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
119790ace30eSBarry Smith   rowners    = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners);
1198ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
119990ace30eSBarry Smith   rowners[0] = 0;
120090ace30eSBarry Smith   for ( i=2; i<=size; i++ ) {
120190ace30eSBarry Smith     rowners[i] += rowners[i-1];
120290ace30eSBarry Smith   }
120390ace30eSBarry Smith 
120490ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
120590ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
120690ace30eSBarry Smith 
120790ace30eSBarry Smith   if (!rank) {
120890ace30eSBarry Smith     vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals);
120990ace30eSBarry Smith 
121090ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
12110752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
121290ace30eSBarry Smith 
121390ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
121490ace30eSBarry Smith     vals_ptr = vals;
121590ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
121690ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
121790ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
121890ace30eSBarry Smith       }
121990ace30eSBarry Smith     }
122090ace30eSBarry Smith 
122190ace30eSBarry Smith     /* read in other processors and ship out */
122290ace30eSBarry Smith     for ( i=1; i<size; i++ ) {
122390ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12240752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1225ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
122690ace30eSBarry Smith     }
12273a40ed3dSBarry Smith   } else {
122890ace30eSBarry Smith     /* receive numeric values */
122990ace30eSBarry Smith     vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals);
123090ace30eSBarry Smith 
123190ace30eSBarry Smith     /* receive message of values*/
1232ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
123390ace30eSBarry Smith 
123490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
123590ace30eSBarry Smith     vals_ptr = vals;
123690ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
123790ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
123890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
123990ace30eSBarry Smith       }
124090ace30eSBarry Smith     }
124190ace30eSBarry Smith   }
1242606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1243606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
12446d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12456d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12463a40ed3dSBarry Smith   PetscFunctionReturn(0);
124790ace30eSBarry Smith }
124890ace30eSBarry Smith 
124990ace30eSBarry Smith 
12505615d1e5SSatish Balay #undef __FUNC__
12515615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense"
125219bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat)
12538965ea79SLois Curfman McInnes {
12548965ea79SLois Curfman McInnes   Mat          A;
12558965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
125619bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
12578965ea79SLois Curfman McInnes   MPI_Status   status;
12588965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
12598965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
126019bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
12613a40ed3dSBarry Smith   int          i, nz, ierr, j,rstart, rend, fd;
12628965ea79SLois Curfman McInnes 
12633a40ed3dSBarry Smith   PetscFunctionBegin;
1264d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1265d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12668965ea79SLois Curfman McInnes   if (!rank) {
126790ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
12680752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1269a8c6a408SBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object");
12708965ea79SLois Curfman McInnes   }
12718965ea79SLois Curfman McInnes 
1272ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
127390ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
127490ace30eSBarry Smith 
127590ace30eSBarry Smith   /*
127690ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
127790ace30eSBarry Smith   */
127890ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
12793a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
12803a40ed3dSBarry Smith     PetscFunctionReturn(0);
128190ace30eSBarry Smith   }
128290ace30eSBarry Smith 
12838965ea79SLois Curfman McInnes   /* determine ownership of all rows */
12848965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
12850452661fSBarry Smith   rowners    = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners);
1286ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
12878965ea79SLois Curfman McInnes   rowners[0] = 0;
12888965ea79SLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
12898965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
12908965ea79SLois Curfman McInnes   }
12918965ea79SLois Curfman McInnes   rstart = rowners[rank];
12928965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
12938965ea79SLois Curfman McInnes 
12948965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
12950452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) );CHKPTRQ(ourlens);
12968965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
12978965ea79SLois Curfman McInnes   if (!rank) {
12980452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) );CHKPTRQ(rowlengths);
12990752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
13000452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(sndcounts);
13018965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
1302ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1303606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1304ca161407SBarry Smith   } else {
1305ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr);
13068965ea79SLois Curfman McInnes   }
13078965ea79SLois Curfman McInnes 
13088965ea79SLois Curfman McInnes   if (!rank) {
13098965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
13100452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(procsnz);
1311549d3d68SSatish Balay     ierr    = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
13128965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
13138965ea79SLois Curfman McInnes       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
13148965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
13158965ea79SLois Curfman McInnes       }
13168965ea79SLois Curfman McInnes     }
1317606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
13188965ea79SLois Curfman McInnes 
13198965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
13208965ea79SLois Curfman McInnes     maxnz = 0;
13218965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
13220452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
13238965ea79SLois Curfman McInnes     }
13240452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) );CHKPTRQ(cols);
13258965ea79SLois Curfman McInnes 
13268965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13278965ea79SLois Curfman McInnes     nz = procsnz[0];
13280452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols);
13290752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13308965ea79SLois Curfman McInnes 
13318965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13328965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
13338965ea79SLois Curfman McInnes       nz   = procsnz[i];
13340752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1335ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13368965ea79SLois Curfman McInnes     }
1337606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13383a40ed3dSBarry Smith   } else {
13398965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13408965ea79SLois Curfman McInnes     nz = 0;
13418965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
13428965ea79SLois Curfman McInnes       nz += ourlens[i];
13438965ea79SLois Curfman McInnes     }
13440452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols);
13458965ea79SLois Curfman McInnes 
13468965ea79SLois Curfman McInnes     /* receive message of column indices*/
1347ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1348ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
1349a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
13508965ea79SLois Curfman McInnes   }
13518965ea79SLois Curfman McInnes 
13528965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1353549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
13548965ea79SLois Curfman McInnes   jj = 0;
13558965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
13568965ea79SLois Curfman McInnes     for ( j=0; j<ourlens[i]; j++ ) {
13578965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
13588965ea79SLois Curfman McInnes       jj++;
13598965ea79SLois Curfman McInnes     }
13608965ea79SLois Curfman McInnes   }
13618965ea79SLois Curfman McInnes 
13628965ea79SLois Curfman McInnes   /* create our matrix */
13638965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
13648965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
13658965ea79SLois Curfman McInnes   }
1366b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
13678965ea79SLois Curfman McInnes   A = *newmat;
13688965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
13698965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
13708965ea79SLois Curfman McInnes   }
13718965ea79SLois Curfman McInnes 
13728965ea79SLois Curfman McInnes   if (!rank) {
13730452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) );CHKPTRQ(vals);
13748965ea79SLois Curfman McInnes 
13758965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
13768965ea79SLois Curfman McInnes     nz = procsnz[0];
13770752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
13788965ea79SLois Curfman McInnes 
13798965ea79SLois Curfman McInnes     /* insert into matrix */
13808965ea79SLois Curfman McInnes     jj      = rstart;
13818965ea79SLois Curfman McInnes     smycols = mycols;
13828965ea79SLois Curfman McInnes     svals   = vals;
13838965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
13848965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13858965ea79SLois Curfman McInnes       smycols += ourlens[i];
13868965ea79SLois Curfman McInnes       svals   += ourlens[i];
13878965ea79SLois Curfman McInnes       jj++;
13888965ea79SLois Curfman McInnes     }
13898965ea79SLois Curfman McInnes 
13908965ea79SLois Curfman McInnes     /* read in other processors and ship out */
13918965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
13928965ea79SLois Curfman McInnes       nz   = procsnz[i];
13930752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1394ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
13958965ea79SLois Curfman McInnes     }
1396606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
13973a40ed3dSBarry Smith   } else {
13988965ea79SLois Curfman McInnes     /* receive numeric values */
13990452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) );CHKPTRQ(vals);
14008965ea79SLois Curfman McInnes 
14018965ea79SLois Curfman McInnes     /* receive message of values*/
1402ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1403ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
1404a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
14058965ea79SLois Curfman McInnes 
14068965ea79SLois Curfman McInnes     /* insert into matrix */
14078965ea79SLois Curfman McInnes     jj      = rstart;
14088965ea79SLois Curfman McInnes     smycols = mycols;
14098965ea79SLois Curfman McInnes     svals   = vals;
14108965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
14118965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14128965ea79SLois Curfman McInnes       smycols += ourlens[i];
14138965ea79SLois Curfman McInnes       svals   += ourlens[i];
14148965ea79SLois Curfman McInnes       jj++;
14158965ea79SLois Curfman McInnes     }
14168965ea79SLois Curfman McInnes   }
1417606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1418606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1419606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1420606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
14218965ea79SLois Curfman McInnes 
14226d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14236d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14243a40ed3dSBarry Smith   PetscFunctionReturn(0);
14258965ea79SLois Curfman McInnes }
142690ace30eSBarry Smith 
142790ace30eSBarry Smith 
142890ace30eSBarry Smith 
142990ace30eSBarry Smith 
143090ace30eSBarry Smith 
1431