xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 35d8aa7f675fa4ee5c80aa031e65d3ef307a4891)
1*35d8aa7fSBarry Smith /*$Id: mpidense.c,v 1.145 2000/10/24 20:25:30 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;
16273d9f13SBarry Smith   int          m = A->m,rstart = mdn->rstart,rank,ierr;
170de54da6SSatish Balay   Scalar       *array;
180de54da6SSatish Balay   MPI_Comm     comm;
190de54da6SSatish Balay 
200de54da6SSatish Balay   PetscFunctionBegin;
21273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
220de54da6SSatish Balay 
230de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
240de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
250de54da6SSatish Balay 
260de54da6SSatish Balay   ierr = 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;
47273d9f13SBarry 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;
52273d9f13SBarry 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;
60273d9f13SBarry 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");
87273d9f13SBarry 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");
92273d9f13SBarry 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;
272*35d8aa7fSBarry Smith   int            *procs,*nprocs,j,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;
280*35d8aa7fSBarry Smith   PetscTruth     found;
2818965ea79SLois Curfman McInnes 
2823a40ed3dSBarry Smith   PetscFunctionBegin;
283b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
2848965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
2858965ea79SLois Curfman McInnes 
2868965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
2870452661fSBarry Smith   nprocs = (int*)PetscMalloc(2*size*sizeof(int));CHKPTRQ(nprocs);
288549d3d68SSatish Balay   ierr   = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
289549d3d68SSatish Balay   procs  = nprocs + size;
2900452661fSBarry Smith   owner  = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(owner); /* see note*/
2918965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
2928965ea79SLois Curfman McInnes     idx = rows[i];
293*35d8aa7fSBarry Smith     found = PETSC_FALSE;
2948965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
2958965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
296*35d8aa7fSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; found = PETSC_TRUE; break;
2978965ea79SLois Curfman McInnes       }
2988965ea79SLois Curfman McInnes     }
29929bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3008965ea79SLois Curfman McInnes   }
3018965ea79SLois Curfman McInnes   nsends = 0;  for (i=0; i<size; i++) { nsends += procs[i];}
3028965ea79SLois Curfman McInnes 
3038965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
3046831982aSBarry Smith   work   = (int*)PetscMalloc(2*size*sizeof(int));CHKPTRQ(work);
3056831982aSBarry Smith   ierr   = MPI_Allreduce(nprocs,work,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr);
3068965ea79SLois Curfman McInnes   nmax   = work[rank];
3076831982aSBarry Smith   nrecvs = work[size+rank];
308606d414cSSatish Balay   ierr   = PetscFree(work);CHKERRQ(ierr);
3098965ea79SLois Curfman McInnes 
3108965ea79SLois Curfman McInnes   /* post receives:   */
3113a40ed3dSBarry Smith   rvalues    = (int*)PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues);
3123a40ed3dSBarry Smith   recv_waits = (MPI_Request*)PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits);
3138965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
314ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3158965ea79SLois Curfman McInnes   }
3168965ea79SLois Curfman McInnes 
3178965ea79SLois Curfman McInnes   /* do sends:
3188965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3198965ea79SLois Curfman McInnes          the ith processor
3208965ea79SLois Curfman McInnes   */
3210452661fSBarry Smith   svalues    = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(svalues);
3227056b6fcSBarry Smith   send_waits = (MPI_Request*)PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
3230452661fSBarry Smith   starts     = (int*)PetscMalloc((size+1)*sizeof(int));CHKPTRQ(starts);
3248965ea79SLois Curfman McInnes   starts[0]  = 0;
3258965ea79SLois Curfman McInnes   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];}
3268965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3278965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3288965ea79SLois Curfman McInnes   }
3298965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3308965ea79SLois Curfman McInnes 
3318965ea79SLois Curfman McInnes   starts[0] = 0;
3328965ea79SLois Curfman McInnes   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[i-1];}
3338965ea79SLois Curfman McInnes   count = 0;
3348965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
3358965ea79SLois Curfman McInnes     if (procs[i]) {
336ca161407SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3378965ea79SLois Curfman McInnes     }
3388965ea79SLois Curfman McInnes   }
339606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3408965ea79SLois Curfman McInnes 
3418965ea79SLois Curfman McInnes   base = owners[rank];
3428965ea79SLois Curfman McInnes 
3438965ea79SLois Curfman McInnes   /*  wait on receives */
3440452661fSBarry Smith   lens   = (int*)PetscMalloc(2*(nrecvs+1)*sizeof(int));CHKPTRQ(lens);
3458965ea79SLois Curfman McInnes   source = lens + nrecvs;
3468965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3478965ea79SLois Curfman McInnes   while (count) {
348ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3498965ea79SLois Curfman McInnes     /* unpack receives into our local space */
350ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3518965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3528965ea79SLois Curfman McInnes     lens[imdex]  = n;
3538965ea79SLois Curfman McInnes     slen += n;
3548965ea79SLois Curfman McInnes     count--;
3558965ea79SLois Curfman McInnes   }
356606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3578965ea79SLois Curfman McInnes 
3588965ea79SLois Curfman McInnes   /* move the data into the send scatter */
3590452661fSBarry Smith   lrows = (int*)PetscMalloc((slen+1)*sizeof(int));CHKPTRQ(lrows);
3608965ea79SLois Curfman McInnes   count = 0;
3618965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
3628965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3638965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
3648965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3658965ea79SLois Curfman McInnes     }
3668965ea79SLois Curfman McInnes   }
367606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
368606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
369606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
370606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3718965ea79SLois Curfman McInnes 
3728965ea79SLois Curfman McInnes   /* actually zap the local rows */
373029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
3748965ea79SLois Curfman McInnes   PLogObjectParent(A,istmp);
375606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3768965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3788965ea79SLois Curfman McInnes 
3798965ea79SLois Curfman McInnes   /* wait on sends */
3808965ea79SLois Curfman McInnes   if (nsends) {
3817056b6fcSBarry Smith     send_status = (MPI_Status*)PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
382ca161407SBarry Smith     ierr        = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
383606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
3848965ea79SLois Curfman McInnes   }
385606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
386606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
3878965ea79SLois Curfman McInnes 
3883a40ed3dSBarry Smith   PetscFunctionReturn(0);
3898965ea79SLois Curfman McInnes }
3908965ea79SLois Curfman McInnes 
3915615d1e5SSatish Balay #undef __FUNC__
392b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMult_MPIDense"
3938f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3948965ea79SLois Curfman McInnes {
39539ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
3968965ea79SLois Curfman McInnes   int          ierr;
397c456f294SBarry Smith 
3983a40ed3dSBarry Smith   PetscFunctionBegin;
39943a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40043a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40144cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4023a40ed3dSBarry Smith   PetscFunctionReturn(0);
4038965ea79SLois Curfman McInnes }
4048965ea79SLois Curfman McInnes 
4055615d1e5SSatish Balay #undef __FUNC__
406b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultAdd_MPIDense"
4078f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4088965ea79SLois Curfman McInnes {
40939ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4108965ea79SLois Curfman McInnes   int          ierr;
411c456f294SBarry Smith 
4123a40ed3dSBarry Smith   PetscFunctionBegin;
41343a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41443a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41544cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4163a40ed3dSBarry Smith   PetscFunctionReturn(0);
4178965ea79SLois Curfman McInnes }
4188965ea79SLois Curfman McInnes 
4195615d1e5SSatish Balay #undef __FUNC__
420b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultTranspose_MPIDense"
4217c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
422096963f5SLois Curfman McInnes {
423096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
424096963f5SLois Curfman McInnes   int          ierr;
4253501a2bdSLois Curfman McInnes   Scalar       zero = 0.0;
426096963f5SLois Curfman McInnes 
4273a40ed3dSBarry Smith   PetscFunctionBegin;
4283501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
4297c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
430537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
431537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4323a40ed3dSBarry Smith   PetscFunctionReturn(0);
433096963f5SLois Curfman McInnes }
434096963f5SLois Curfman McInnes 
4355615d1e5SSatish Balay #undef __FUNC__
436b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultTransposeAdd_MPIDense"
4377c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
438096963f5SLois Curfman McInnes {
439096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
440096963f5SLois Curfman McInnes   int          ierr;
441096963f5SLois Curfman McInnes 
4423a40ed3dSBarry Smith   PetscFunctionBegin;
4433501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4447c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
445537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
446537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4473a40ed3dSBarry Smith   PetscFunctionReturn(0);
448096963f5SLois Curfman McInnes }
449096963f5SLois Curfman McInnes 
4505615d1e5SSatish Balay #undef __FUNC__
451b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetDiagonal_MPIDense"
4528f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v)
4538965ea79SLois Curfman McInnes {
45439ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
455096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
456273d9f13SBarry Smith   int          ierr,len,i,n,m = A->m,radd;
45744cd7ae7SLois Curfman McInnes   Scalar       *x,zero = 0.0;
458ed3cc1f0SBarry Smith 
4593a40ed3dSBarry Smith   PetscFunctionBegin;
460273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
461096963f5SLois Curfman McInnes   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
462096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
463273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
464273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
4657ddc982cSLois Curfman McInnes   radd = a->rstart*m;
46644cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
467096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
468096963f5SLois Curfman McInnes   }
4699a8c540fSSatish Balay   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
4703a40ed3dSBarry Smith   PetscFunctionReturn(0);
4718965ea79SLois Curfman McInnes }
4728965ea79SLois Curfman McInnes 
4735615d1e5SSatish Balay #undef __FUNC__
474b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDestroy_MPIDense"
475e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat)
4768965ea79SLois Curfman McInnes {
4773501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4788965ea79SLois Curfman McInnes   int          ierr;
479ed3cc1f0SBarry Smith 
4803a40ed3dSBarry Smith   PetscFunctionBegin;
48194d884c6SBarry Smith 
482aa482453SBarry Smith #if defined(PETSC_USE_LOG)
483273d9f13SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N);
4848965ea79SLois Curfman McInnes #endif
4858798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
486606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
4873501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
4883501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4893501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
490622d7880SLois Curfman McInnes   if (mdn->factor) {
491606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
492606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
493606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
494606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
495622d7880SLois Curfman McInnes   }
496606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
4973a40ed3dSBarry Smith   PetscFunctionReturn(0);
4988965ea79SLois Curfman McInnes }
49939ddd567SLois Curfman McInnes 
5005615d1e5SSatish Balay #undef __FUNC__
501b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense_Binary"
50239ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer)
5038965ea79SLois Curfman McInnes {
50439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
5058965ea79SLois Curfman McInnes   int          ierr;
5067056b6fcSBarry Smith 
5073a40ed3dSBarry Smith   PetscFunctionBegin;
50839ddd567SLois Curfman McInnes   if (mdn->size == 1) {
50939ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5108965ea79SLois Curfman McInnes   }
51129bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5123a40ed3dSBarry Smith   PetscFunctionReturn(0);
5138965ea79SLois Curfman McInnes }
5148965ea79SLois Curfman McInnes 
5155615d1e5SSatish Balay #undef __FUNC__
516b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense_ASCIIorDraworSocket"
517f1af5d2fSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,Viewer viewer)
5188965ea79SLois Curfman McInnes {
51939ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
52077ed5343SBarry Smith   int          ierr,format,size = mdn->size,rank = mdn->rank;
52119bcc07fSBarry Smith   ViewerType   vtype;
522f1af5d2fSBarry Smith   PetscTruth   isascii,isdraw;
523b9b97703SBarry Smith   Viewer       sviewer;
5248965ea79SLois Curfman McInnes 
5253a40ed3dSBarry Smith   PetscFunctionBegin;
526f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
527f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
528f1af5d2fSBarry Smith   if (isascii) {
5293a40ed3dSBarry Smith     ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
530888f2ed8SSatish Balay     ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
531639f9d9dSBarry Smith     if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
5324e220ebcSLois Curfman McInnes       MatInfo info;
533888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
534273d9f13SBarry Smith       ierr = ViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m,
5356831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
5366831982aSBarry Smith       ierr = ViewerFlush(viewer);CHKERRQ(ierr);
5373501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5383a40ed3dSBarry Smith       PetscFunctionReturn(0);
53996f6c058SBarry Smith     } else if (format == VIEWER_FORMAT_ASCII_INFO) {
5403a40ed3dSBarry Smith       PetscFunctionReturn(0);
5418965ea79SLois Curfman McInnes     }
542f1af5d2fSBarry Smith   } else if (isdraw) {
543f1af5d2fSBarry Smith     Draw       draw;
544f1af5d2fSBarry Smith     PetscTruth isnull;
545f1af5d2fSBarry Smith 
546f1af5d2fSBarry Smith     ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
547f1af5d2fSBarry Smith     ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr);
548f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
549f1af5d2fSBarry Smith   }
55077ed5343SBarry Smith 
5518965ea79SLois Curfman McInnes   if (size == 1) {
55239ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5533a40ed3dSBarry Smith   } else {
5548965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5558965ea79SLois Curfman McInnes     Mat          A;
556273d9f13SBarry Smith     int          M = mat->M,N = mat->N,m,row,i,nz,*cols;
55739ddd567SLois Curfman McInnes     Scalar       *vals;
5588965ea79SLois Curfman McInnes 
5598965ea79SLois Curfman McInnes     if (!rank) {
560f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5613a40ed3dSBarry Smith     } else {
562f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5638965ea79SLois Curfman McInnes     }
5648965ea79SLois Curfman McInnes     PLogObjectParent(mat,A);
5658965ea79SLois Curfman McInnes 
56639ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
56739ddd567SLois Curfman McInnes        but it's quick for now */
568273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
5698965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
57039ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57139ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
57239ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57339ddd567SLois Curfman McInnes       row++;
5748965ea79SLois Curfman McInnes     }
5758965ea79SLois Curfman McInnes 
5766d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5776d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5786831982aSBarry Smith     ierr = ViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
579b9b97703SBarry Smith     if (!rank) {
5806831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5818965ea79SLois Curfman McInnes     }
582b9b97703SBarry Smith     ierr = ViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
5836831982aSBarry Smith     ierr = ViewerFlush(viewer);CHKERRQ(ierr);
5848965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
5858965ea79SLois Curfman McInnes   }
5863a40ed3dSBarry Smith   PetscFunctionReturn(0);
5878965ea79SLois Curfman McInnes }
5888965ea79SLois Curfman McInnes 
5895615d1e5SSatish Balay #undef __FUNC__
590b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense"
591e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer)
5928965ea79SLois Curfman McInnes {
59339ddd567SLois Curfman McInnes   int        ierr;
594f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
5958965ea79SLois Curfman McInnes 
596433994e6SBarry Smith   PetscFunctionBegin;
5970f5bd95cSBarry Smith 
5986831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
5996831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
600f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
601f1af5d2fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
6020f5bd95cSBarry Smith 
603f1af5d2fSBarry Smith   if (isascii || issocket || isdraw) {
604f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6050f5bd95cSBarry Smith   } else if (isbinary) {
6063a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6075cd90555SBarry Smith   } else {
60829bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6098965ea79SLois Curfman McInnes   }
6103a40ed3dSBarry Smith   PetscFunctionReturn(0);
6118965ea79SLois Curfman McInnes }
6128965ea79SLois Curfman McInnes 
6135615d1e5SSatish Balay #undef __FUNC__
614b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetInfo_MPIDense"
6158f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6168965ea79SLois Curfman McInnes {
6173501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6183501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6194e220ebcSLois Curfman McInnes   int          ierr;
620329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6218965ea79SLois Curfman McInnes 
6223a40ed3dSBarry Smith   PetscFunctionBegin;
623273d9f13SBarry Smith   info->rows_global    = (double)A->M;
624273d9f13SBarry Smith   info->columns_global = (double)A->N;
625273d9f13SBarry Smith   info->rows_local     = (double)A->m;
626273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6274e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6284e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6294e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6304e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6318965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6324e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6334e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6344e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6354e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6364e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6378965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
638f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
6394e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6404e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6414e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6424e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6434e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6448965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
645f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
6464e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6474e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6484e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6494e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6504e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6518965ea79SLois Curfman McInnes   }
6524e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6534e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6544e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6553a40ed3dSBarry Smith   PetscFunctionReturn(0);
6568965ea79SLois Curfman McInnes }
6578965ea79SLois Curfman McInnes 
6585615d1e5SSatish Balay #undef __FUNC__
659b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_MPIDense"
6608f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6618965ea79SLois Curfman McInnes {
66239ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
663273d9f13SBarry Smith   int          ierr;
6648965ea79SLois Curfman McInnes 
6653a40ed3dSBarry Smith   PetscFunctionBegin;
6666d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6676d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
6684787f768SSatish Balay       op == MAT_NEW_NONZERO_LOCATION_ERR ||
6694787f768SSatish Balay       op == MAT_NEW_NONZERO_ALLOCATION_ERR ||
670219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
671219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED) {
672273d9f13SBarry Smith         ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
673b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
674273d9f13SBarry Smith         a->roworiented = PETSC_TRUE;
675273d9f13SBarry Smith         ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
676b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
677219d9a1aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
6786d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
6796d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
680b51ba29fSSatish Balay              op == MAT_YES_NEW_DIAGONALS ||
681b51ba29fSSatish Balay              op == MAT_USE_HASH_TABLE) {
682981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
6833a40ed3dSBarry Smith   } else if (op == MAT_COLUMN_ORIENTED) {
684273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
685273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
6863782ba37SSatish Balay   } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) {
687273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
6883a40ed3dSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS) {
68929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
6903a40ed3dSBarry Smith   } else {
69129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
6923a40ed3dSBarry Smith   }
6933a40ed3dSBarry Smith   PetscFunctionReturn(0);
6948965ea79SLois Curfman McInnes }
6958965ea79SLois Curfman McInnes 
6965615d1e5SSatish Balay #undef __FUNC__
697b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_MPIDense"
6988f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
6998965ea79SLois Curfman McInnes {
7003501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7013a40ed3dSBarry Smith 
7023a40ed3dSBarry Smith   PetscFunctionBegin;
7034c49b128SBarry Smith   if (m) *m = mat->rstart;
7044c49b128SBarry Smith   if (n) *n = mat->rend;
7053a40ed3dSBarry Smith   PetscFunctionReturn(0);
7068965ea79SLois Curfman McInnes }
7078965ea79SLois Curfman McInnes 
7085615d1e5SSatish Balay #undef __FUNC__
709b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_MPIDense"
7108f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
7118965ea79SLois Curfman McInnes {
7123501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7133a40ed3dSBarry Smith   int          lrow,rstart = mat->rstart,rend = mat->rend,ierr;
7148965ea79SLois Curfman McInnes 
7153a40ed3dSBarry Smith   PetscFunctionBegin;
71629bbc08cSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
7178965ea79SLois Curfman McInnes   lrow = row - rstart;
7183a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7193a40ed3dSBarry Smith   PetscFunctionReturn(0);
7208965ea79SLois Curfman McInnes }
7218965ea79SLois Curfman McInnes 
7225615d1e5SSatish Balay #undef __FUNC__
723b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_MPIDense"
7248f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
7258965ea79SLois Curfman McInnes {
726606d414cSSatish Balay   int ierr;
727606d414cSSatish Balay 
7283a40ed3dSBarry Smith   PetscFunctionBegin;
729606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
730606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7313a40ed3dSBarry Smith   PetscFunctionReturn(0);
7328965ea79SLois Curfman McInnes }
7338965ea79SLois Curfman McInnes 
7345615d1e5SSatish Balay #undef __FUNC__
735b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDiagonalScale_MPIDense"
7365b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7375b2fa520SLois Curfman McInnes {
7385b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7395b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7405b2fa520SLois Curfman McInnes   Scalar       *l,*r,x,*v;
741273d9f13SBarry Smith   int          ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7425b2fa520SLois Curfman McInnes 
7435b2fa520SLois Curfman McInnes   PetscFunctionBegin;
74472d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7455b2fa520SLois Curfman McInnes   if (ll) {
74672d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
74729bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
7485b2fa520SLois Curfman McInnes     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7495b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7505b2fa520SLois Curfman McInnes       x = l[i];
7515b2fa520SLois Curfman McInnes       v = mat->v + i;
7525b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7535b2fa520SLois Curfman McInnes     }
7545b2fa520SLois Curfman McInnes     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
7555b2fa520SLois Curfman McInnes     PLogFlops(n*m);
7565b2fa520SLois Curfman McInnes   }
7575b2fa520SLois Curfman McInnes   if (rr) {
75872d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
75929bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7605b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7615b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7625b2fa520SLois Curfman McInnes     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
7635b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7645b2fa520SLois Curfman McInnes       x = r[i];
7655b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7665b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7675b2fa520SLois Curfman McInnes     }
76872d926a5SLois Curfman McInnes     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
7695b2fa520SLois Curfman McInnes     PLogFlops(n*m);
7705b2fa520SLois Curfman McInnes   }
7715b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7725b2fa520SLois Curfman McInnes }
7735b2fa520SLois Curfman McInnes 
7745b2fa520SLois Curfman McInnes #undef __FUNC__
775b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatNorm_MPIDense"
776329f5518SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *norm)
777096963f5SLois Curfman McInnes {
7783501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7793501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7803501a2bdSLois Curfman McInnes   int          ierr,i,j;
781329f5518SBarry Smith   PetscReal    sum = 0.0;
7823501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
7833501a2bdSLois Curfman McInnes 
7843a40ed3dSBarry Smith   PetscFunctionBegin;
7853501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
7863501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm);CHKERRQ(ierr);
7873501a2bdSLois Curfman McInnes   } else {
7883501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
789273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
790aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
791329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
7923501a2bdSLois Curfman McInnes #else
7933501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7943501a2bdSLois Curfman McInnes #endif
7953501a2bdSLois Curfman McInnes       }
796ca161407SBarry Smith       ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
7973501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
798273d9f13SBarry Smith       PLogFlops(2*mdn->A->n*mdn->A->m);
7993a40ed3dSBarry Smith     } else if (type == NORM_1) {
800329f5518SBarry Smith       PetscReal *tmp,*tmp2;
801273d9f13SBarry Smith       tmp  = (PetscReal*)PetscMalloc(2*A->N*sizeof(PetscReal));CHKPTRQ(tmp);
802273d9f13SBarry Smith       tmp2 = tmp + A->N;
803273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
804096963f5SLois Curfman McInnes       *norm = 0.0;
8053501a2bdSLois Curfman McInnes       v = mat->v;
806273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
807273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
80867e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8093501a2bdSLois Curfman McInnes         }
8103501a2bdSLois Curfman McInnes       }
811273d9f13SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
812273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
8133501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
8143501a2bdSLois Curfman McInnes       }
815606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
816273d9f13SBarry Smith       PLogFlops(A->n*A->m);
8173a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
818329f5518SBarry Smith       PetscReal ntemp;
8193501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
820ca161407SBarry Smith       ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
8213a40ed3dSBarry Smith     } else {
82229bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8233501a2bdSLois Curfman McInnes     }
8243501a2bdSLois Curfman McInnes   }
8253a40ed3dSBarry Smith   PetscFunctionReturn(0);
8263501a2bdSLois Curfman McInnes }
8273501a2bdSLois Curfman McInnes 
8285615d1e5SSatish Balay #undef __FUNC__
829b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatTranspose_MPIDense"
8308f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8313501a2bdSLois Curfman McInnes {
8323501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8333501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8343501a2bdSLois Curfman McInnes   Mat          B;
835273d9f13SBarry Smith   int          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8363501a2bdSLois Curfman McInnes   int          j,i,ierr;
8373501a2bdSLois Curfman McInnes   Scalar       *v;
8383501a2bdSLois Curfman McInnes 
8393a40ed3dSBarry Smith   PetscFunctionBegin;
8407c922b88SBarry Smith   if (!matout && M != N) {
84129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8427056b6fcSBarry Smith   }
8437056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8443501a2bdSLois Curfman McInnes 
845273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
8460452661fSBarry Smith   rwork = (int*)PetscMalloc(n*sizeof(int));CHKPTRQ(rwork);
8473501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8483501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8493501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8503501a2bdSLois Curfman McInnes     v   += m;
8513501a2bdSLois Curfman McInnes   }
852606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8536d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8546d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8557c922b88SBarry Smith   if (matout) {
8563501a2bdSLois Curfman McInnes     *matout = B;
8573501a2bdSLois Curfman McInnes   } else {
858273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8593501a2bdSLois Curfman McInnes   }
8603a40ed3dSBarry Smith   PetscFunctionReturn(0);
861096963f5SLois Curfman McInnes }
862096963f5SLois Curfman McInnes 
863d9eff348SSatish Balay #include "petscblaslapack.h"
8645615d1e5SSatish Balay #undef __FUNC__
865b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatScale_MPIDense"
8668f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA)
86744cd7ae7SLois Curfman McInnes {
86844cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
86944cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
87044cd7ae7SLois Curfman McInnes   int          one = 1,nz;
87144cd7ae7SLois Curfman McInnes 
8723a40ed3dSBarry Smith   PetscFunctionBegin;
873273d9f13SBarry Smith   nz = inA->m*inA->N;
87444cd7ae7SLois Curfman McInnes   BLscal_(&nz,alpha,a->v,&one);
87544cd7ae7SLois Curfman McInnes   PLogFlops(nz);
8763a40ed3dSBarry Smith   PetscFunctionReturn(0);
87744cd7ae7SLois Curfman McInnes }
87844cd7ae7SLois Curfman McInnes 
8795609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
880ca44d042SBarry Smith EXTERN int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **);
8818965ea79SLois Curfman McInnes 
882273d9f13SBarry Smith #undef __FUNC__
883273d9f13SBarry Smith #define __FUNC__ /*<a name="MatSetUpPreallocation_MPIDense"></a>*/"MatSetUpPreallocation_MPIDense"
884273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A)
885273d9f13SBarry Smith {
886273d9f13SBarry Smith   int        ierr;
887273d9f13SBarry Smith 
888273d9f13SBarry Smith   PetscFunctionBegin;
889273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
890273d9f13SBarry Smith   PetscFunctionReturn(0);
891273d9f13SBarry Smith }
892273d9f13SBarry Smith 
8938965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
89409dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
89509dc0095SBarry Smith        MatGetRow_MPIDense,
89609dc0095SBarry Smith        MatRestoreRow_MPIDense,
89709dc0095SBarry Smith        MatMult_MPIDense,
89809dc0095SBarry Smith        MatMultAdd_MPIDense,
8997c922b88SBarry Smith        MatMultTranspose_MPIDense,
9007c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9018965ea79SLois Curfman McInnes        0,
90209dc0095SBarry Smith        0,
90309dc0095SBarry Smith        0,
90409dc0095SBarry Smith        0,
90509dc0095SBarry Smith        0,
90609dc0095SBarry Smith        0,
90709dc0095SBarry Smith        0,
90809dc0095SBarry Smith        MatTranspose_MPIDense,
90909dc0095SBarry Smith        MatGetInfo_MPIDense,0,
91009dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9115b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
91209dc0095SBarry Smith        MatNorm_MPIDense,
91309dc0095SBarry Smith        MatAssemblyBegin_MPIDense,
91409dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
91509dc0095SBarry Smith        0,
91609dc0095SBarry Smith        MatSetOption_MPIDense,
91709dc0095SBarry Smith        MatZeroEntries_MPIDense,
91809dc0095SBarry Smith        MatZeroRows_MPIDense,
91909dc0095SBarry Smith        0,
92009dc0095SBarry Smith        0,
92109dc0095SBarry Smith        0,
92209dc0095SBarry Smith        0,
923273d9f13SBarry Smith        MatSetUpPreallocation_MPIDense,
924273d9f13SBarry Smith        0,
92539ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
92609dc0095SBarry Smith        0,
92709dc0095SBarry Smith        0,
92809dc0095SBarry Smith        MatGetArray_MPIDense,
92909dc0095SBarry Smith        MatRestoreArray_MPIDense,
9305609ef8eSBarry Smith        MatDuplicate_MPIDense,
93109dc0095SBarry Smith        0,
93209dc0095SBarry Smith        0,
93309dc0095SBarry Smith        0,
93409dc0095SBarry Smith        0,
93509dc0095SBarry Smith        0,
9362ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
93709dc0095SBarry Smith        0,
93809dc0095SBarry Smith        MatGetValues_MPIDense,
93909dc0095SBarry Smith        0,
94009dc0095SBarry Smith        0,
94109dc0095SBarry Smith        MatScale_MPIDense,
94209dc0095SBarry Smith        0,
94309dc0095SBarry Smith        0,
94409dc0095SBarry Smith        0,
94509dc0095SBarry Smith        MatGetBlockSize_MPIDense,
94609dc0095SBarry Smith        0,
94709dc0095SBarry Smith        0,
94809dc0095SBarry Smith        0,
94909dc0095SBarry Smith        0,
95009dc0095SBarry Smith        0,
95109dc0095SBarry Smith        0,
95209dc0095SBarry Smith        0,
95309dc0095SBarry Smith        0,
95409dc0095SBarry Smith        0,
955ca3fa75bSLois Curfman McInnes        MatGetSubMatrix_MPIDense,
956b9b97703SBarry Smith        MatDestroy_MPIDense,
957b9b97703SBarry Smith        MatView_MPIDense,
95809dc0095SBarry Smith        MatGetMaps_Petsc};
9598965ea79SLois Curfman McInnes 
960273d9f13SBarry Smith EXTERN_C_BEGIN
961273d9f13SBarry Smith #undef __FUNC__
962273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreate_MPIDense"
963273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat)
964273d9f13SBarry Smith {
965273d9f13SBarry Smith   Mat_MPIDense *a;
966273d9f13SBarry Smith   int          ierr,i;
967273d9f13SBarry Smith 
968273d9f13SBarry Smith   PetscFunctionBegin;
969273d9f13SBarry Smith   mat->data         = (void*)(a = PetscNew(Mat_MPIDense));CHKPTRQ(a);
970273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
971273d9f13SBarry Smith   mat->factor       = 0;
972273d9f13SBarry Smith   mat->mapping      = 0;
973273d9f13SBarry Smith 
974273d9f13SBarry Smith   a->factor       = 0;
975273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
976273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
977273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
978273d9f13SBarry Smith 
979273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
980273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
981273d9f13SBarry Smith   a->nvec = mat->n;
982273d9f13SBarry Smith 
983273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
984273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
985273d9f13SBarry Smith   ierr = MapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
986273d9f13SBarry Smith   ierr = MapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
987273d9f13SBarry Smith 
988273d9f13SBarry Smith   /* build local table of row and column ownerships */
989273d9f13SBarry Smith   a->rowners = (int*)PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners);
990273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
991273d9f13SBarry Smith   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
992273d9f13SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
993273d9f13SBarry Smith   a->rowners[0] = 0;
994273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
995273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
996273d9f13SBarry Smith   }
997273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
998273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
999273d9f13SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1000273d9f13SBarry Smith   a->cowners[0] = 0;
1001273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1002273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1003273d9f13SBarry Smith   }
1004273d9f13SBarry Smith 
1005273d9f13SBarry Smith   /* build cache for off array entries formed */
1006273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1007273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1008273d9f13SBarry Smith 
1009273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1010273d9f13SBarry Smith   a->lvec        = 0;
1011273d9f13SBarry Smith   a->Mvctx       = 0;
1012273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1013273d9f13SBarry Smith 
1014273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1015273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1016273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1017273d9f13SBarry Smith   PetscFunctionReturn(0);
1018273d9f13SBarry Smith }
1019273d9f13SBarry Smith EXTERN_C_END
1020273d9f13SBarry Smith 
1021273d9f13SBarry Smith #undef __FUNC__
1022273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMPIDenseSetPreallocation"
1023273d9f13SBarry Smith /*@C
1024273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1025273d9f13SBarry Smith 
1026273d9f13SBarry Smith    Not collective
1027273d9f13SBarry Smith 
1028273d9f13SBarry Smith    Input Parameters:
1029273d9f13SBarry Smith .  A - the matrix
1030273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1031273d9f13SBarry Smith    to control all matrix memory allocation.
1032273d9f13SBarry Smith 
1033273d9f13SBarry Smith    Notes:
1034273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1035273d9f13SBarry Smith    storage by columns.
1036273d9f13SBarry Smith 
1037273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1038273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1039273d9f13SBarry Smith    set data=PETSC_NULL.
1040273d9f13SBarry Smith 
1041273d9f13SBarry Smith    Level: intermediate
1042273d9f13SBarry Smith 
1043273d9f13SBarry Smith .keywords: matrix,dense, parallel
1044273d9f13SBarry Smith 
1045273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1046273d9f13SBarry Smith @*/
1047273d9f13SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,Scalar *data)
1048273d9f13SBarry Smith {
1049273d9f13SBarry Smith   Mat_MPIDense *a;
1050273d9f13SBarry Smith   int          ierr;
1051273d9f13SBarry Smith   PetscTruth   flg2;
1052273d9f13SBarry Smith 
1053273d9f13SBarry Smith   PetscFunctionBegin;
1054273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)mat,MATMPIDENSE,&flg2);CHKERRQ(ierr);
1055273d9f13SBarry Smith   if (!flg2) PetscFunctionReturn(0);
1056273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
1057273d9f13SBarry Smith   /* Note:  For now, when data is specified above, this assumes the user correctly
1058273d9f13SBarry Smith    allocates the local dense storage space.  We should add error checking. */
1059273d9f13SBarry Smith 
1060273d9f13SBarry Smith   a    = (Mat_MPIDense*)mat->data;
1061273d9f13SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mat->m,mat->N,data,&a->A);CHKERRQ(ierr);
1062273d9f13SBarry Smith   PLogObjectParent(mat,a->A);
1063273d9f13SBarry Smith   PetscFunctionReturn(0);
1064273d9f13SBarry Smith }
1065273d9f13SBarry Smith 
10665615d1e5SSatish Balay #undef __FUNC__
1067b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateMPIDense"
10688965ea79SLois Curfman McInnes /*@C
106939ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
10708965ea79SLois Curfman McInnes 
1071db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1072db81eaa0SLois Curfman McInnes 
10738965ea79SLois Curfman McInnes    Input Parameters:
1074db81eaa0SLois Curfman McInnes +  comm - MPI communicator
10758965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1076db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
10778965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1078db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1079db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1080dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
10818965ea79SLois Curfman McInnes 
10828965ea79SLois Curfman McInnes    Output Parameter:
1083477f1c0bSLois Curfman McInnes .  A - the matrix
10848965ea79SLois Curfman McInnes 
1085b259b22eSLois Curfman McInnes    Notes:
108639ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
108739ddd567SLois Curfman McInnes    storage by columns.
10888965ea79SLois Curfman McInnes 
108918f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
109018f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1091b4fd4287SBarry Smith    set data=PETSC_NULL.
109218f449edSLois Curfman McInnes 
10938965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
10948965ea79SLois Curfman McInnes    (possibly both).
10958965ea79SLois Curfman McInnes 
1096027ccd11SLois Curfman McInnes    Level: intermediate
1097027ccd11SLois Curfman McInnes 
109839ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
10998965ea79SLois Curfman McInnes 
110039ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
11018965ea79SLois Curfman McInnes @*/
1102477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
11038965ea79SLois Curfman McInnes {
1104273d9f13SBarry Smith   int ierr,size;
11058965ea79SLois Curfman McInnes 
11063a40ed3dSBarry Smith   PetscFunctionBegin;
1107273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1108273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1109273d9f13SBarry Smith   if (size > 1) {
1110273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1111273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1112273d9f13SBarry Smith   } else {
1113273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1114273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
11158c469469SLois Curfman McInnes   }
11163a40ed3dSBarry Smith   PetscFunctionReturn(0);
11178965ea79SLois Curfman McInnes }
11188965ea79SLois Curfman McInnes 
11195615d1e5SSatish Balay #undef __FUNC__
1120b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDuplicate_MPIDense"
11215609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11228965ea79SLois Curfman McInnes {
11238965ea79SLois Curfman McInnes   Mat          mat;
11243501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
112539ddd567SLois Curfman McInnes   int          ierr;
11262ba99913SLois Curfman McInnes   FactorCtx    *factor;
11278965ea79SLois Curfman McInnes 
11283a40ed3dSBarry Smith   PetscFunctionBegin;
11298965ea79SLois Curfman McInnes   *newmat       = 0;
1130273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1131273d9f13SBarry Smith   ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr);
11320452661fSBarry Smith   mat->data         = (void*)(a = PetscNew(Mat_MPIDense));CHKPTRQ(a);
1133549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
11343501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1135c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1136273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
11378965ea79SLois Curfman McInnes 
11382ba99913SLois Curfman McInnes   if (oldmat->factor) {
11392ba99913SLois Curfman McInnes     a->factor = (FactorCtx*)(factor = PetscNew(FactorCtx));CHKPTRQ(factor);
11402ba99913SLois Curfman McInnes     /* copy factor contents ... add this code! */
11412ba99913SLois Curfman McInnes   } else a->factor = 0;
11428965ea79SLois Curfman McInnes 
11438965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
11448965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
11458965ea79SLois Curfman McInnes   a->size         = oldmat->size;
11468965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1147e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1148b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
11493782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
11500452661fSBarry Smith   a->rowners = (int*)PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners);
1151f09e8eb9SSatish Balay   PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1152549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
11538798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
11548965ea79SLois Curfman McInnes 
1155329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
11565609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
11578965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
11588965ea79SLois Curfman McInnes   *newmat = mat;
11593a40ed3dSBarry Smith   PetscFunctionReturn(0);
11608965ea79SLois Curfman McInnes }
11618965ea79SLois Curfman McInnes 
1162e090d566SSatish Balay #include "petscsys.h"
11638965ea79SLois Curfman McInnes 
11645615d1e5SSatish Balay #undef __FUNC__
1165b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_MPIDense_DenseInFile"
116690ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat)
116790ace30eSBarry Smith {
116840011551SBarry Smith   int        *rowners,i,size,rank,m,ierr,nz,j;
116990ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
117090ace30eSBarry Smith   MPI_Status status;
117190ace30eSBarry Smith 
11723a40ed3dSBarry Smith   PetscFunctionBegin;
1173d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1174d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
117590ace30eSBarry Smith 
117690ace30eSBarry Smith   /* determine ownership of all rows */
117790ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
117890ace30eSBarry Smith   rowners    = (int*)PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners);
1179ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
118090ace30eSBarry Smith   rowners[0] = 0;
118190ace30eSBarry Smith   for (i=2; i<=size; i++) {
118290ace30eSBarry Smith     rowners[i] += rowners[i-1];
118390ace30eSBarry Smith   }
118490ace30eSBarry Smith 
118590ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
118690ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
118790ace30eSBarry Smith 
118890ace30eSBarry Smith   if (!rank) {
118990ace30eSBarry Smith     vals = (Scalar*)PetscMalloc(m*N*sizeof(Scalar));CHKPTRQ(vals);
119090ace30eSBarry Smith 
119190ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
11920752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
119390ace30eSBarry Smith 
119490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
119590ace30eSBarry Smith     vals_ptr = vals;
119690ace30eSBarry Smith     for (i=0; i<m; i++) {
119790ace30eSBarry Smith       for (j=0; j<N; j++) {
119890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
119990ace30eSBarry Smith       }
120090ace30eSBarry Smith     }
120190ace30eSBarry Smith 
120290ace30eSBarry Smith     /* read in other processors and ship out */
120390ace30eSBarry Smith     for (i=1; i<size; i++) {
120490ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12050752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1206ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
120790ace30eSBarry Smith     }
12083a40ed3dSBarry Smith   } else {
120990ace30eSBarry Smith     /* receive numeric values */
121090ace30eSBarry Smith     vals = (Scalar*)PetscMalloc(m*N*sizeof(Scalar));CHKPTRQ(vals);
121190ace30eSBarry Smith 
121290ace30eSBarry Smith     /* receive message of values*/
1213ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
121490ace30eSBarry Smith 
121590ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
121690ace30eSBarry Smith     vals_ptr = vals;
121790ace30eSBarry Smith     for (i=0; i<m; i++) {
121890ace30eSBarry Smith       for (j=0; j<N; j++) {
121990ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
122090ace30eSBarry Smith       }
122190ace30eSBarry Smith     }
122290ace30eSBarry Smith   }
1223606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1224606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
12256d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12266d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12273a40ed3dSBarry Smith   PetscFunctionReturn(0);
122890ace30eSBarry Smith }
122990ace30eSBarry Smith 
1230273d9f13SBarry Smith EXTERN_C_BEGIN
12315615d1e5SSatish Balay #undef __FUNC__
1232b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_MPIDense"
123319bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat)
12348965ea79SLois Curfman McInnes {
12358965ea79SLois Curfman McInnes   Mat          A;
12368965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
123719bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
12388965ea79SLois Curfman McInnes   MPI_Status   status;
12398965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
12408965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols;
124119bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
12423a40ed3dSBarry Smith   int          i,nz,ierr,j,rstart,rend,fd;
12438965ea79SLois Curfman McInnes 
12443a40ed3dSBarry Smith   PetscFunctionBegin;
1245d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1246d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12478965ea79SLois Curfman McInnes   if (!rank) {
124890ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
12490752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
125029bbc08cSBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
12518965ea79SLois Curfman McInnes   }
12528965ea79SLois Curfman McInnes 
1253ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
125490ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
125590ace30eSBarry Smith 
125690ace30eSBarry Smith   /*
125790ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
125890ace30eSBarry Smith   */
125990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
12603a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
12613a40ed3dSBarry Smith     PetscFunctionReturn(0);
126290ace30eSBarry Smith   }
126390ace30eSBarry Smith 
12648965ea79SLois Curfman McInnes   /* determine ownership of all rows */
12658965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
12660452661fSBarry Smith   rowners    = (int*)PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners);
1267ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
12688965ea79SLois Curfman McInnes   rowners[0] = 0;
12698965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
12708965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
12718965ea79SLois Curfman McInnes   }
12728965ea79SLois Curfman McInnes   rstart = rowners[rank];
12738965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
12748965ea79SLois Curfman McInnes 
12758965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
12760452661fSBarry Smith   ourlens = (int*)PetscMalloc(2*(rend-rstart)*sizeof(int));CHKPTRQ(ourlens);
12778965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
12788965ea79SLois Curfman McInnes   if (!rank) {
12790452661fSBarry Smith     rowlengths = (int*)PetscMalloc(M*sizeof(int));CHKPTRQ(rowlengths);
12800752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
12810452661fSBarry Smith     sndcounts = (int*)PetscMalloc(size*sizeof(int));CHKPTRQ(sndcounts);
12828965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1283ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1284606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1285ca161407SBarry Smith   } else {
1286ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
12878965ea79SLois Curfman McInnes   }
12888965ea79SLois Curfman McInnes 
12898965ea79SLois Curfman McInnes   if (!rank) {
12908965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
12910452661fSBarry Smith     procsnz = (int*)PetscMalloc(size*sizeof(int));CHKPTRQ(procsnz);
1292549d3d68SSatish Balay     ierr    = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
12938965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
12948965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
12958965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
12968965ea79SLois Curfman McInnes       }
12978965ea79SLois Curfman McInnes     }
1298606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
12998965ea79SLois Curfman McInnes 
13008965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
13018965ea79SLois Curfman McInnes     maxnz = 0;
13028965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13030452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
13048965ea79SLois Curfman McInnes     }
13050452661fSBarry Smith     cols = (int*)PetscMalloc(maxnz*sizeof(int));CHKPTRQ(cols);
13068965ea79SLois Curfman McInnes 
13078965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13088965ea79SLois Curfman McInnes     nz = procsnz[0];
13090452661fSBarry Smith     mycols = (int*)PetscMalloc(nz*sizeof(int));CHKPTRQ(mycols);
13100752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13118965ea79SLois Curfman McInnes 
13128965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13138965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13148965ea79SLois Curfman McInnes       nz   = procsnz[i];
13150752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1316ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13178965ea79SLois Curfman McInnes     }
1318606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13193a40ed3dSBarry Smith   } else {
13208965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13218965ea79SLois Curfman McInnes     nz = 0;
13228965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13238965ea79SLois Curfman McInnes       nz += ourlens[i];
13248965ea79SLois Curfman McInnes     }
13250452661fSBarry Smith     mycols = (int*)PetscMalloc(nz*sizeof(int));CHKPTRQ(mycols);
13268965ea79SLois Curfman McInnes 
13278965ea79SLois Curfman McInnes     /* receive message of column indices*/
1328ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1329ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
133029bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13318965ea79SLois Curfman McInnes   }
13328965ea79SLois Curfman McInnes 
13338965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1334549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
13358965ea79SLois Curfman McInnes   jj = 0;
13368965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13378965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
13388965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
13398965ea79SLois Curfman McInnes       jj++;
13408965ea79SLois Curfman McInnes     }
13418965ea79SLois Curfman McInnes   }
13428965ea79SLois Curfman McInnes 
13438965ea79SLois Curfman McInnes   /* create our matrix */
13448965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13458965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
13468965ea79SLois Curfman McInnes   }
1347b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
13488965ea79SLois Curfman McInnes   A = *newmat;
13498965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13508965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
13518965ea79SLois Curfman McInnes   }
13528965ea79SLois Curfman McInnes 
13538965ea79SLois Curfman McInnes   if (!rank) {
13540452661fSBarry Smith     vals = (Scalar*)PetscMalloc(maxnz*sizeof(Scalar));CHKPTRQ(vals);
13558965ea79SLois Curfman McInnes 
13568965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
13578965ea79SLois Curfman McInnes     nz = procsnz[0];
13580752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
13598965ea79SLois Curfman McInnes 
13608965ea79SLois Curfman McInnes     /* insert into matrix */
13618965ea79SLois Curfman McInnes     jj      = rstart;
13628965ea79SLois Curfman McInnes     smycols = mycols;
13638965ea79SLois Curfman McInnes     svals   = vals;
13648965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13658965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13668965ea79SLois Curfman McInnes       smycols += ourlens[i];
13678965ea79SLois Curfman McInnes       svals   += ourlens[i];
13688965ea79SLois Curfman McInnes       jj++;
13698965ea79SLois Curfman McInnes     }
13708965ea79SLois Curfman McInnes 
13718965ea79SLois Curfman McInnes     /* read in other processors and ship out */
13728965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13738965ea79SLois Curfman McInnes       nz   = procsnz[i];
13740752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1375ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
13768965ea79SLois Curfman McInnes     }
1377606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
13783a40ed3dSBarry Smith   } else {
13798965ea79SLois Curfman McInnes     /* receive numeric values */
13800452661fSBarry Smith     vals = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(vals);
13818965ea79SLois Curfman McInnes 
13828965ea79SLois Curfman McInnes     /* receive message of values*/
1383ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1384ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
138529bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13868965ea79SLois Curfman McInnes 
13878965ea79SLois Curfman McInnes     /* insert into matrix */
13888965ea79SLois Curfman McInnes     jj      = rstart;
13898965ea79SLois Curfman McInnes     smycols = mycols;
13908965ea79SLois Curfman McInnes     svals   = vals;
13918965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13928965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13938965ea79SLois Curfman McInnes       smycols += ourlens[i];
13948965ea79SLois Curfman McInnes       svals   += ourlens[i];
13958965ea79SLois Curfman McInnes       jj++;
13968965ea79SLois Curfman McInnes     }
13978965ea79SLois Curfman McInnes   }
1398606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1399606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1400606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1401606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
14028965ea79SLois Curfman McInnes 
14036d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14046d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14053a40ed3dSBarry Smith   PetscFunctionReturn(0);
14068965ea79SLois Curfman McInnes }
1407273d9f13SBarry Smith EXTERN_C_END
140890ace30eSBarry Smith 
140990ace30eSBarry Smith 
141090ace30eSBarry Smith 
141190ace30eSBarry Smith 
141290ace30eSBarry Smith 
1413