xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 273d9f13de75c4ed17021f7f2c11eebb99d26f0d)
1*273d9f13SBarry Smith /*$Id: mpidense.c,v 1.144 2000/09/28 21:10:58 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__
12b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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*273d9f13SBarry Smith   int          m = A->m,rstart = mdn->rstart,rank,ierr;
170de54da6SSatish Balay   Scalar       *array;
180de54da6SSatish Balay   MPI_Comm     comm;
190de54da6SSatish Balay 
200de54da6SSatish Balay   PetscFunctionBegin;
21*273d9f13SBarry 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 = 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 
39ca44d042SBarry Smith EXTERN int MatSetUpMultiply_MPIDense(Mat);
407ef1d9bdSSatish Balay 
415615d1e5SSatish Balay #undef __FUNC__
42b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
47*273d9f13SBarry Smith   PetscTruth   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;
52*273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"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;
60*273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"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__
78b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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++) {
8629bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
87*273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"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++) {
9129bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
92*273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
9329bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"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 {
9829bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"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__
105b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
117b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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);
130b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
131b9b97703SBarry Smith   ierr = ISGetLocalSize(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) {
14229bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"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 */
14732828cfdSBarry Smith     ierr = MatCreateMPIDense(A->comm,nrows,cs,PETSC_DECIDE,ncols,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__
173b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
181b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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)) {
19329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"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);
199c38d4ed2SBarry 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__
204b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
240b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
252b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
267b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
282b9b97703SBarry Smith   ierr = ISGetLocalSize(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     }
29829bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"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__
391b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
405b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
419b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultTranspose_MPIDense"
4207c922b88SBarry Smith int MatMultTranspose_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);
4287c922b88SBarry Smith   ierr = MatMultTranspose_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__
435b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultTransposeAdd_MPIDense"
4367c922b88SBarry Smith int MatMultTransposeAdd_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);
4437c922b88SBarry Smith   ierr = MatMultTranspose_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__
450b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
455*273d9f13SBarry Smith   int          ierr,len,i,n,m = A->m,radd;
45644cd7ae7SLois Curfman McInnes   Scalar       *x,zero = 0.0;
457ed3cc1f0SBarry Smith 
4583a40ed3dSBarry Smith   PetscFunctionBegin;
459*273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
460096963f5SLois Curfman McInnes   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
461096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
462*273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
463*273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->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__
473b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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 
481aa482453SBarry Smith #if defined(PETSC_USE_LOG)
482*273d9f13SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N);
4838965ea79SLois Curfman McInnes #endif
4848798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
485606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
4863501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
4873501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4883501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
489622d7880SLois Curfman McInnes   if (mdn->factor) {
490606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
491606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
492606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
493606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
494622d7880SLois Curfman McInnes   }
495606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
4963a40ed3dSBarry Smith   PetscFunctionReturn(0);
4978965ea79SLois Curfman McInnes }
49839ddd567SLois Curfman McInnes 
4995615d1e5SSatish Balay #undef __FUNC__
500b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense_Binary"
50139ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer)
5028965ea79SLois Curfman McInnes {
50339ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
5048965ea79SLois Curfman McInnes   int          ierr;
5057056b6fcSBarry Smith 
5063a40ed3dSBarry Smith   PetscFunctionBegin;
50739ddd567SLois Curfman McInnes   if (mdn->size == 1) {
50839ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5098965ea79SLois Curfman McInnes   }
51029bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5113a40ed3dSBarry Smith   PetscFunctionReturn(0);
5128965ea79SLois Curfman McInnes }
5138965ea79SLois Curfman McInnes 
5145615d1e5SSatish Balay #undef __FUNC__
515b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense_ASCIIorDraworSocket"
516f1af5d2fSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,Viewer viewer)
5178965ea79SLois Curfman McInnes {
51839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
51977ed5343SBarry Smith   int          ierr,format,size = mdn->size,rank = mdn->rank;
52019bcc07fSBarry Smith   ViewerType   vtype;
521f1af5d2fSBarry Smith   PetscTruth   isascii,isdraw;
522b9b97703SBarry Smith   Viewer       sviewer;
5238965ea79SLois Curfman McInnes 
5243a40ed3dSBarry Smith   PetscFunctionBegin;
525f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
526f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
527f1af5d2fSBarry Smith   if (isascii) {
5283a40ed3dSBarry Smith     ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
529888f2ed8SSatish Balay     ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
530639f9d9dSBarry Smith     if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
5314e220ebcSLois Curfman McInnes       MatInfo info;
532888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
533*273d9f13SBarry Smith       ierr = ViewerASCIISynchronizedPrintf(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);
5356831982aSBarry Smith       ierr = ViewerFlush(viewer);CHKERRQ(ierr);
5363501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5373a40ed3dSBarry Smith       PetscFunctionReturn(0);
53896f6c058SBarry Smith     } else if (format == VIEWER_FORMAT_ASCII_INFO) {
5393a40ed3dSBarry Smith       PetscFunctionReturn(0);
5408965ea79SLois Curfman McInnes     }
541f1af5d2fSBarry Smith   } else if (isdraw) {
542f1af5d2fSBarry Smith     Draw       draw;
543f1af5d2fSBarry Smith     PetscTruth isnull;
544f1af5d2fSBarry Smith 
545f1af5d2fSBarry Smith     ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
546f1af5d2fSBarry Smith     ierr = DrawIsNull(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;
555*273d9f13SBarry 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     }
5638965ea79SLois Curfman McInnes     PLogObjectParent(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 */
567*273d9f13SBarry 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);
5776831982aSBarry Smith     ierr = ViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
578b9b97703SBarry Smith     if (!rank) {
5796831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5808965ea79SLois Curfman McInnes     }
581b9b97703SBarry Smith     ierr = ViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
5826831982aSBarry Smith     ierr = ViewerFlush(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__
589b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense"
590e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer)
5918965ea79SLois Curfman McInnes {
59239ddd567SLois Curfman McInnes   int        ierr;
593f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
5948965ea79SLois Curfman McInnes 
595433994e6SBarry Smith   PetscFunctionBegin;
5960f5bd95cSBarry Smith 
5976831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
5986831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
599f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
600f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&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__
613b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
622*273d9f13SBarry Smith   info->rows_global    = (double)A->M;
623*273d9f13SBarry Smith   info->columns_global = (double)A->N;
624*273d9f13SBarry Smith   info->rows_local     = (double)A->m;
625*273d9f13SBarry 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__
658b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
662*273d9f13SBarry 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) {
671*273d9f13SBarry Smith         ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
672b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
673*273d9f13SBarry Smith         a->roworiented = PETSC_TRUE;
674*273d9f13SBarry 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) {
681981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
6823a40ed3dSBarry Smith   } else if (op == MAT_COLUMN_ORIENTED) {
683*273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
684*273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
6853782ba37SSatish Balay   } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) {
686*273d9f13SBarry 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__
696b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
708b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
722b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
734b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
740*273d9f13SBarry 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);
7545b2fa520SLois Curfman McInnes     PLogFlops(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);
7685b2fa520SLois Curfman McInnes     PLogFlops(n*m);
7695b2fa520SLois Curfman McInnes   }
7705b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7715b2fa520SLois Curfman McInnes }
7725b2fa520SLois Curfman McInnes 
7735b2fa520SLois Curfman McInnes #undef __FUNC__
774b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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) {
788*273d9f13SBarry 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);
797*273d9f13SBarry Smith       PLogFlops(2*mdn->A->n*mdn->A->m);
7983a40ed3dSBarry Smith     } else if (type == NORM_1) {
799329f5518SBarry Smith       PetscReal *tmp,*tmp2;
800*273d9f13SBarry Smith       tmp  = (PetscReal*)PetscMalloc(2*A->N*sizeof(PetscReal));CHKPTRQ(tmp);
801*273d9f13SBarry Smith       tmp2 = tmp + A->N;
802*273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
803096963f5SLois Curfman McInnes       *norm = 0.0;
8043501a2bdSLois Curfman McInnes       v = mat->v;
805*273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
806*273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
80767e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8083501a2bdSLois Curfman McInnes         }
8093501a2bdSLois Curfman McInnes       }
810*273d9f13SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
811*273d9f13SBarry 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);
815*273d9f13SBarry Smith       PLogFlops(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__
828b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
834*273d9f13SBarry 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 
844*273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
8450452661fSBarry Smith   rwork = (int*)PetscMalloc(n*sizeof(int));CHKPTRQ(rwork);
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 {
857*273d9f13SBarry 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__
864b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
872*273d9f13SBarry Smith   nz = inA->m*inA->N;
87344cd7ae7SLois Curfman McInnes   BLscal_(&nz,alpha,a->v,&one);
87444cd7ae7SLois Curfman McInnes   PLogFlops(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 
881*273d9f13SBarry Smith #undef __FUNC__
882*273d9f13SBarry Smith #define __FUNC__ /*<a name="MatSetUpPreallocation_MPIDense"></a>*/"MatSetUpPreallocation_MPIDense"
883*273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A)
884*273d9f13SBarry Smith {
885*273d9f13SBarry Smith   int        ierr;
886*273d9f13SBarry Smith 
887*273d9f13SBarry Smith   PetscFunctionBegin;
888*273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
889*273d9f13SBarry Smith   PetscFunctionReturn(0);
890*273d9f13SBarry Smith }
891*273d9f13SBarry 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,
922*273d9f13SBarry Smith        MatSetUpPreallocation_MPIDense,
923*273d9f13SBarry 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 
959*273d9f13SBarry Smith EXTERN_C_BEGIN
960*273d9f13SBarry Smith #undef __FUNC__
961*273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreate_MPIDense"
962*273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat)
963*273d9f13SBarry Smith {
964*273d9f13SBarry Smith   Mat_MPIDense *a;
965*273d9f13SBarry Smith   int          ierr,i;
966*273d9f13SBarry Smith 
967*273d9f13SBarry Smith   PetscFunctionBegin;
968*273d9f13SBarry Smith   mat->data         = (void*)(a = PetscNew(Mat_MPIDense));CHKPTRQ(a);
969*273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
970*273d9f13SBarry Smith   mat->factor       = 0;
971*273d9f13SBarry Smith   mat->mapping      = 0;
972*273d9f13SBarry Smith 
973*273d9f13SBarry Smith   a->factor       = 0;
974*273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
975*273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
976*273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
977*273d9f13SBarry Smith 
978*273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
979*273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
980*273d9f13SBarry Smith   a->nvec = mat->n;
981*273d9f13SBarry Smith 
982*273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
983*273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
984*273d9f13SBarry Smith   ierr = MapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
985*273d9f13SBarry Smith   ierr = MapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
986*273d9f13SBarry Smith 
987*273d9f13SBarry Smith   /* build local table of row and column ownerships */
988*273d9f13SBarry Smith   a->rowners = (int*)PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners);
989*273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
990*273d9f13SBarry Smith   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
991*273d9f13SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
992*273d9f13SBarry Smith   a->rowners[0] = 0;
993*273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
994*273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
995*273d9f13SBarry Smith   }
996*273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
997*273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
998*273d9f13SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
999*273d9f13SBarry Smith   a->cowners[0] = 0;
1000*273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1001*273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1002*273d9f13SBarry Smith   }
1003*273d9f13SBarry Smith 
1004*273d9f13SBarry Smith   /* build cache for off array entries formed */
1005*273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1006*273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1007*273d9f13SBarry Smith 
1008*273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1009*273d9f13SBarry Smith   a->lvec        = 0;
1010*273d9f13SBarry Smith   a->Mvctx       = 0;
1011*273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1012*273d9f13SBarry Smith 
1013*273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1014*273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1015*273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1016*273d9f13SBarry Smith   PetscFunctionReturn(0);
1017*273d9f13SBarry Smith }
1018*273d9f13SBarry Smith EXTERN_C_END
1019*273d9f13SBarry Smith 
1020*273d9f13SBarry Smith #undef __FUNC__
1021*273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMPIDenseSetPreallocation"
1022*273d9f13SBarry Smith /*@C
1023*273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1024*273d9f13SBarry Smith 
1025*273d9f13SBarry Smith    Not collective
1026*273d9f13SBarry Smith 
1027*273d9f13SBarry Smith    Input Parameters:
1028*273d9f13SBarry Smith .  A - the matrix
1029*273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1030*273d9f13SBarry Smith    to control all matrix memory allocation.
1031*273d9f13SBarry Smith 
1032*273d9f13SBarry Smith    Notes:
1033*273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1034*273d9f13SBarry Smith    storage by columns.
1035*273d9f13SBarry Smith 
1036*273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1037*273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1038*273d9f13SBarry Smith    set data=PETSC_NULL.
1039*273d9f13SBarry Smith 
1040*273d9f13SBarry Smith    Level: intermediate
1041*273d9f13SBarry Smith 
1042*273d9f13SBarry Smith .keywords: matrix,dense, parallel
1043*273d9f13SBarry Smith 
1044*273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1045*273d9f13SBarry Smith @*/
1046*273d9f13SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,Scalar *data)
1047*273d9f13SBarry Smith {
1048*273d9f13SBarry Smith   Mat_MPIDense *a;
1049*273d9f13SBarry Smith   int          ierr;
1050*273d9f13SBarry Smith   PetscTruth   flg2;
1051*273d9f13SBarry Smith 
1052*273d9f13SBarry Smith   PetscFunctionBegin;
1053*273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)mat,MATMPIDENSE,&flg2);CHKERRQ(ierr);
1054*273d9f13SBarry Smith   if (!flg2) PetscFunctionReturn(0);
1055*273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
1056*273d9f13SBarry Smith   /* Note:  For now, when data is specified above, this assumes the user correctly
1057*273d9f13SBarry Smith    allocates the local dense storage space.  We should add error checking. */
1058*273d9f13SBarry Smith 
1059*273d9f13SBarry Smith   a    = (Mat_MPIDense*)mat->data;
1060*273d9f13SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mat->m,mat->N,data,&a->A);CHKERRQ(ierr);
1061*273d9f13SBarry Smith   PLogObjectParent(mat,a->A);
1062*273d9f13SBarry Smith   PetscFunctionReturn(0);
1063*273d9f13SBarry Smith }
1064*273d9f13SBarry Smith 
10655615d1e5SSatish Balay #undef __FUNC__
1066b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateMPIDense"
10678965ea79SLois Curfman McInnes /*@C
106839ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
10698965ea79SLois Curfman McInnes 
1070db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1071db81eaa0SLois Curfman McInnes 
10728965ea79SLois Curfman McInnes    Input Parameters:
1073db81eaa0SLois Curfman McInnes +  comm - MPI communicator
10748965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1075db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
10768965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1077db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1078db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1079dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
10808965ea79SLois Curfman McInnes 
10818965ea79SLois Curfman McInnes    Output Parameter:
1082477f1c0bSLois Curfman McInnes .  A - the matrix
10838965ea79SLois Curfman McInnes 
1084b259b22eSLois Curfman McInnes    Notes:
108539ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
108639ddd567SLois Curfman McInnes    storage by columns.
10878965ea79SLois Curfman McInnes 
108818f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
108918f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1090b4fd4287SBarry Smith    set data=PETSC_NULL.
109118f449edSLois Curfman McInnes 
10928965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
10938965ea79SLois Curfman McInnes    (possibly both).
10948965ea79SLois Curfman McInnes 
1095027ccd11SLois Curfman McInnes    Level: intermediate
1096027ccd11SLois Curfman McInnes 
109739ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
10988965ea79SLois Curfman McInnes 
109939ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
11008965ea79SLois Curfman McInnes @*/
1101477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
11028965ea79SLois Curfman McInnes {
1103*273d9f13SBarry Smith   int ierr,size;
11048965ea79SLois Curfman McInnes 
11053a40ed3dSBarry Smith   PetscFunctionBegin;
1106*273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1107*273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1108*273d9f13SBarry Smith   if (size > 1) {
1109*273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1110*273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1111*273d9f13SBarry Smith   } else {
1112*273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1113*273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
11148c469469SLois Curfman McInnes   }
11153a40ed3dSBarry Smith   PetscFunctionReturn(0);
11168965ea79SLois Curfman McInnes }
11178965ea79SLois Curfman McInnes 
11185615d1e5SSatish Balay #undef __FUNC__
1119b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDuplicate_MPIDense"
11205609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11218965ea79SLois Curfman McInnes {
11228965ea79SLois Curfman McInnes   Mat          mat;
11233501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
112439ddd567SLois Curfman McInnes   int          ierr;
11252ba99913SLois Curfman McInnes   FactorCtx    *factor;
11268965ea79SLois Curfman McInnes 
11273a40ed3dSBarry Smith   PetscFunctionBegin;
11288965ea79SLois Curfman McInnes   *newmat       = 0;
1129*273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1130*273d9f13SBarry Smith   ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr);
11310452661fSBarry Smith   mat->data         = (void*)(a = PetscNew(Mat_MPIDense));CHKPTRQ(a);
1132549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
11333501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1134c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1135*273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
11368965ea79SLois Curfman McInnes 
11372ba99913SLois Curfman McInnes   if (oldmat->factor) {
11382ba99913SLois Curfman McInnes     a->factor = (FactorCtx*)(factor = PetscNew(FactorCtx));CHKPTRQ(factor);
11392ba99913SLois Curfman McInnes     /* copy factor contents ... add this code! */
11402ba99913SLois Curfman McInnes   } else a->factor = 0;
11418965ea79SLois Curfman McInnes 
11428965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
11438965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
11448965ea79SLois Curfman McInnes   a->size         = oldmat->size;
11458965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1146e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1147b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
11483782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
11490452661fSBarry Smith   a->rowners = (int*)PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners);
1150f09e8eb9SSatish Balay   PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1151549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
11528798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
11538965ea79SLois Curfman McInnes 
1154329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
11555609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
11568965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
11578965ea79SLois Curfman McInnes   *newmat = mat;
11583a40ed3dSBarry Smith   PetscFunctionReturn(0);
11598965ea79SLois Curfman McInnes }
11608965ea79SLois Curfman McInnes 
1161e090d566SSatish Balay #include "petscsys.h"
11628965ea79SLois Curfman McInnes 
11635615d1e5SSatish Balay #undef __FUNC__
1164b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_MPIDense_DenseInFile"
116590ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat)
116690ace30eSBarry Smith {
116740011551SBarry Smith   int        *rowners,i,size,rank,m,ierr,nz,j;
116890ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
116990ace30eSBarry Smith   MPI_Status status;
117090ace30eSBarry Smith 
11713a40ed3dSBarry Smith   PetscFunctionBegin;
1172d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1173d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
117490ace30eSBarry Smith 
117590ace30eSBarry Smith   /* determine ownership of all rows */
117690ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
117790ace30eSBarry Smith   rowners    = (int*)PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners);
1178ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
117990ace30eSBarry Smith   rowners[0] = 0;
118090ace30eSBarry Smith   for (i=2; i<=size; i++) {
118190ace30eSBarry Smith     rowners[i] += rowners[i-1];
118290ace30eSBarry Smith   }
118390ace30eSBarry Smith 
118490ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
118590ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
118690ace30eSBarry Smith 
118790ace30eSBarry Smith   if (!rank) {
118890ace30eSBarry Smith     vals = (Scalar*)PetscMalloc(m*N*sizeof(Scalar));CHKPTRQ(vals);
118990ace30eSBarry Smith 
119090ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
11910752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
119290ace30eSBarry Smith 
119390ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
119490ace30eSBarry Smith     vals_ptr = vals;
119590ace30eSBarry Smith     for (i=0; i<m; i++) {
119690ace30eSBarry Smith       for (j=0; j<N; j++) {
119790ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
119890ace30eSBarry Smith       }
119990ace30eSBarry Smith     }
120090ace30eSBarry Smith 
120190ace30eSBarry Smith     /* read in other processors and ship out */
120290ace30eSBarry Smith     for (i=1; i<size; i++) {
120390ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12040752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1205ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
120690ace30eSBarry Smith     }
12073a40ed3dSBarry Smith   } else {
120890ace30eSBarry Smith     /* receive numeric values */
120990ace30eSBarry Smith     vals = (Scalar*)PetscMalloc(m*N*sizeof(Scalar));CHKPTRQ(vals);
121090ace30eSBarry Smith 
121190ace30eSBarry Smith     /* receive message of values*/
1212ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
121390ace30eSBarry Smith 
121490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
121590ace30eSBarry Smith     vals_ptr = vals;
121690ace30eSBarry Smith     for (i=0; i<m; i++) {
121790ace30eSBarry Smith       for (j=0; j<N; j++) {
121890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
121990ace30eSBarry Smith       }
122090ace30eSBarry Smith     }
122190ace30eSBarry Smith   }
1222606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1223606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
12246d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12256d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12263a40ed3dSBarry Smith   PetscFunctionReturn(0);
122790ace30eSBarry Smith }
122890ace30eSBarry Smith 
1229*273d9f13SBarry Smith EXTERN_C_BEGIN
12305615d1e5SSatish Balay #undef __FUNC__
1231b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_MPIDense"
123219bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat)
12338965ea79SLois Curfman McInnes {
12348965ea79SLois Curfman McInnes   Mat          A;
12358965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
123619bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
12378965ea79SLois Curfman McInnes   MPI_Status   status;
12388965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
12398965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols;
124019bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
12413a40ed3dSBarry Smith   int          i,nz,ierr,j,rstart,rend,fd;
12428965ea79SLois Curfman McInnes 
12433a40ed3dSBarry Smith   PetscFunctionBegin;
1244d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1245d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12468965ea79SLois Curfman McInnes   if (!rank) {
124790ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
12480752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
124929bbc08cSBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
12508965ea79SLois Curfman McInnes   }
12518965ea79SLois Curfman McInnes 
1252ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
125390ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
125490ace30eSBarry Smith 
125590ace30eSBarry Smith   /*
125690ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
125790ace30eSBarry Smith   */
125890ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
12593a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
12603a40ed3dSBarry Smith     PetscFunctionReturn(0);
126190ace30eSBarry Smith   }
126290ace30eSBarry Smith 
12638965ea79SLois Curfman McInnes   /* determine ownership of all rows */
12648965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
12650452661fSBarry Smith   rowners    = (int*)PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners);
1266ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
12678965ea79SLois Curfman McInnes   rowners[0] = 0;
12688965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
12698965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
12708965ea79SLois Curfman McInnes   }
12718965ea79SLois Curfman McInnes   rstart = rowners[rank];
12728965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
12738965ea79SLois Curfman McInnes 
12748965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
12750452661fSBarry Smith   ourlens = (int*)PetscMalloc(2*(rend-rstart)*sizeof(int));CHKPTRQ(ourlens);
12768965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
12778965ea79SLois Curfman McInnes   if (!rank) {
12780452661fSBarry Smith     rowlengths = (int*)PetscMalloc(M*sizeof(int));CHKPTRQ(rowlengths);
12790752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
12800452661fSBarry Smith     sndcounts = (int*)PetscMalloc(size*sizeof(int));CHKPTRQ(sndcounts);
12818965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1282ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1283606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1284ca161407SBarry Smith   } else {
1285ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
12868965ea79SLois Curfman McInnes   }
12878965ea79SLois Curfman McInnes 
12888965ea79SLois Curfman McInnes   if (!rank) {
12898965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
12900452661fSBarry Smith     procsnz = (int*)PetscMalloc(size*sizeof(int));CHKPTRQ(procsnz);
1291549d3d68SSatish Balay     ierr    = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
12928965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
12938965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
12948965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
12958965ea79SLois Curfman McInnes       }
12968965ea79SLois Curfman McInnes     }
1297606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
12988965ea79SLois Curfman McInnes 
12998965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
13008965ea79SLois Curfman McInnes     maxnz = 0;
13018965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13020452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
13038965ea79SLois Curfman McInnes     }
13040452661fSBarry Smith     cols = (int*)PetscMalloc(maxnz*sizeof(int));CHKPTRQ(cols);
13058965ea79SLois Curfman McInnes 
13068965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13078965ea79SLois Curfman McInnes     nz = procsnz[0];
13080452661fSBarry Smith     mycols = (int*)PetscMalloc(nz*sizeof(int));CHKPTRQ(mycols);
13090752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13108965ea79SLois Curfman McInnes 
13118965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13128965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13138965ea79SLois Curfman McInnes       nz   = procsnz[i];
13140752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1315ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13168965ea79SLois Curfman McInnes     }
1317606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13183a40ed3dSBarry Smith   } else {
13198965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13208965ea79SLois Curfman McInnes     nz = 0;
13218965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13228965ea79SLois Curfman McInnes       nz += ourlens[i];
13238965ea79SLois Curfman McInnes     }
13240452661fSBarry Smith     mycols = (int*)PetscMalloc(nz*sizeof(int));CHKPTRQ(mycols);
13258965ea79SLois Curfman McInnes 
13268965ea79SLois Curfman McInnes     /* receive message of column indices*/
1327ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1328ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
132929bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13308965ea79SLois Curfman McInnes   }
13318965ea79SLois Curfman McInnes 
13328965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1333549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
13348965ea79SLois Curfman McInnes   jj = 0;
13358965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13368965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
13378965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
13388965ea79SLois Curfman McInnes       jj++;
13398965ea79SLois Curfman McInnes     }
13408965ea79SLois Curfman McInnes   }
13418965ea79SLois Curfman McInnes 
13428965ea79SLois Curfman McInnes   /* create our matrix */
13438965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13448965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
13458965ea79SLois Curfman McInnes   }
1346b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
13478965ea79SLois Curfman McInnes   A = *newmat;
13488965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13498965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
13508965ea79SLois Curfman McInnes   }
13518965ea79SLois Curfman McInnes 
13528965ea79SLois Curfman McInnes   if (!rank) {
13530452661fSBarry Smith     vals = (Scalar*)PetscMalloc(maxnz*sizeof(Scalar));CHKPTRQ(vals);
13548965ea79SLois Curfman McInnes 
13558965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
13568965ea79SLois Curfman McInnes     nz = procsnz[0];
13570752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
13588965ea79SLois Curfman McInnes 
13598965ea79SLois Curfman McInnes     /* insert into matrix */
13608965ea79SLois Curfman McInnes     jj      = rstart;
13618965ea79SLois Curfman McInnes     smycols = mycols;
13628965ea79SLois Curfman McInnes     svals   = vals;
13638965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13648965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13658965ea79SLois Curfman McInnes       smycols += ourlens[i];
13668965ea79SLois Curfman McInnes       svals   += ourlens[i];
13678965ea79SLois Curfman McInnes       jj++;
13688965ea79SLois Curfman McInnes     }
13698965ea79SLois Curfman McInnes 
13708965ea79SLois Curfman McInnes     /* read in other processors and ship out */
13718965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13728965ea79SLois Curfman McInnes       nz   = procsnz[i];
13730752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1374ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
13758965ea79SLois Curfman McInnes     }
1376606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
13773a40ed3dSBarry Smith   } else {
13788965ea79SLois Curfman McInnes     /* receive numeric values */
13790452661fSBarry Smith     vals = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(vals);
13808965ea79SLois Curfman McInnes 
13818965ea79SLois Curfman McInnes     /* receive message of values*/
1382ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1383ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
138429bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13858965ea79SLois Curfman McInnes 
13868965ea79SLois Curfman McInnes     /* insert into matrix */
13878965ea79SLois Curfman McInnes     jj      = rstart;
13888965ea79SLois Curfman McInnes     smycols = mycols;
13898965ea79SLois Curfman McInnes     svals   = vals;
13908965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13918965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13928965ea79SLois Curfman McInnes       smycols += ourlens[i];
13938965ea79SLois Curfman McInnes       svals   += ourlens[i];
13948965ea79SLois Curfman McInnes       jj++;
13958965ea79SLois Curfman McInnes     }
13968965ea79SLois Curfman McInnes   }
1397606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1398606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1399606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1400606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
14018965ea79SLois Curfman McInnes 
14026d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14036d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14043a40ed3dSBarry Smith   PetscFunctionReturn(0);
14058965ea79SLois Curfman McInnes }
1406*273d9f13SBarry Smith EXTERN_C_END
140790ace30eSBarry Smith 
140890ace30eSBarry Smith 
140990ace30eSBarry Smith 
141090ace30eSBarry Smith 
141190ace30eSBarry Smith 
1412