xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision b1d4fb267e72f6b087ed013bbfbbee1418c3f503)
173f4d377SMatthew Knepley /*$Id: mpidense.c,v 1.159 2001/08/10 03:30:41 bsmith Exp $*/
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"
88965ea79SLois Curfman McInnes 
90de54da6SSatish Balay EXTERN_C_BEGIN
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
120de54da6SSatish Balay int MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
130de54da6SSatish Balay {
140de54da6SSatish Balay   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
15cfce73b9SSatish Balay   int          m = A->m,rstart = mdn->rstart,ierr;
1687828ca2SBarry Smith   PetscScalar  *array;
170de54da6SSatish Balay   MPI_Comm     comm;
180de54da6SSatish Balay 
190de54da6SSatish Balay   PetscFunctionBegin;
20273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
210de54da6SSatish Balay 
220de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
230de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
240de54da6SSatish Balay 
250de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
260de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
270de54da6SSatish Balay   ierr = MatCreateSeqDense(comm,m,m,array+m*rstart,B);CHKERRQ(ierr);
280de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
290de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
300de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
310de54da6SSatish Balay 
320de54da6SSatish Balay   *iscopy = PETSC_TRUE;
330de54da6SSatish Balay   PetscFunctionReturn(0);
340de54da6SSatish Balay }
350de54da6SSatish Balay EXTERN_C_END
360de54da6SSatish Balay 
374a2ae208SSatish Balay #undef __FUNCT__
384a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
39f15d580aSBarry Smith int MatSetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],const PetscScalar v[],InsertMode addv)
408965ea79SLois Curfman McInnes {
4139b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
4239b7565bSBarry Smith   int          ierr,i,j,rstart = A->rstart,rend = A->rend,row;
43273d9f13SBarry Smith   PetscTruth   roworiented = A->roworiented;
448965ea79SLois Curfman McInnes 
453a40ed3dSBarry Smith   PetscFunctionBegin;
468965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
475ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
48273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
498965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
508965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
5139b7565bSBarry Smith       if (roworiented) {
5239b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
533a40ed3dSBarry Smith       } else {
548965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
555ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
56273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
5739b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
5839b7565bSBarry Smith         }
598965ea79SLois Curfman McInnes       }
603a40ed3dSBarry Smith     } else {
613782ba37SSatish Balay       if (!A->donotstash) {
6239b7565bSBarry Smith         if (roworiented) {
638798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
64d36fbae8SSatish Balay         } else {
658798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
6639b7565bSBarry Smith         }
67b49de8d1SLois Curfman McInnes       }
68b49de8d1SLois Curfman McInnes     }
693782ba37SSatish Balay   }
703a40ed3dSBarry Smith   PetscFunctionReturn(0);
71b49de8d1SLois Curfman McInnes }
72b49de8d1SLois Curfman McInnes 
734a2ae208SSatish Balay #undef __FUNCT__
744a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
75f15d580aSBarry Smith int MatGetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],PetscScalar v[])
76b49de8d1SLois Curfman McInnes {
77b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
78b49de8d1SLois Curfman McInnes   int          ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row;
79b49de8d1SLois Curfman McInnes 
803a40ed3dSBarry Smith   PetscFunctionBegin;
81b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
8229bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
83273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
84b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
85b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
86b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
8729bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
88273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
8929bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
90a8c6a408SBarry Smith         }
91b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
92b49de8d1SLois Curfman McInnes       }
93a8c6a408SBarry Smith     } else {
9429bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
958965ea79SLois Curfman McInnes     }
968965ea79SLois Curfman McInnes   }
973a40ed3dSBarry Smith   PetscFunctionReturn(0);
988965ea79SLois Curfman McInnes }
998965ea79SLois Curfman McInnes 
1004a2ae208SSatish Balay #undef __FUNCT__
1014a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
1024e7234bfSBarry Smith int MatGetArray_MPIDense(Mat A,PetscScalar *array[])
103ff14e315SSatish Balay {
104ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
105ff14e315SSatish Balay   int          ierr;
106ff14e315SSatish Balay 
1073a40ed3dSBarry Smith   PetscFunctionBegin;
108ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1093a40ed3dSBarry Smith   PetscFunctionReturn(0);
110ff14e315SSatish Balay }
111ff14e315SSatish Balay 
1124a2ae208SSatish Balay #undef __FUNCT__
1134a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
114ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
115ca3fa75bSLois Curfman McInnes {
116ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
117ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
118cfce73b9SSatish Balay   int          i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
11987828ca2SBarry Smith   PetscScalar  *av,*bv,*v = lmat->v;
120ca3fa75bSLois Curfman McInnes   Mat          newmat;
121ca3fa75bSLois Curfman McInnes 
122ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
123ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
124ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
125b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
126b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
127ca3fa75bSLois Curfman McInnes 
128ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1297eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1307eba5e9cSLois Curfman McInnes      original matrix! */
131ca3fa75bSLois Curfman McInnes 
132ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1337eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
134ca3fa75bSLois Curfman McInnes 
135ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
136ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
13729bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
1387eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
139ca3fa75bSLois Curfman McInnes     newmat = *B;
140ca3fa75bSLois Curfman McInnes   } else {
141ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
14232828cfdSBarry Smith     ierr = MatCreateMPIDense(A->comm,nrows,cs,PETSC_DECIDE,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr);
143ca3fa75bSLois Curfman McInnes   }
144ca3fa75bSLois Curfman McInnes 
145ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
146ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
147ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
148ca3fa75bSLois Curfman McInnes 
149ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
150ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
151ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
1527eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
153ca3fa75bSLois Curfman McInnes     }
154ca3fa75bSLois Curfman McInnes   }
155ca3fa75bSLois Curfman McInnes 
156ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
157ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159ca3fa75bSLois Curfman McInnes 
160ca3fa75bSLois Curfman McInnes   /* Free work space */
161ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
162ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
163ca3fa75bSLois Curfman McInnes   *B = newmat;
164ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
165ca3fa75bSLois Curfman McInnes }
166ca3fa75bSLois Curfman McInnes 
1674a2ae208SSatish Balay #undef __FUNCT__
1684a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
1694e7234bfSBarry Smith int MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
170ff14e315SSatish Balay {
1713a40ed3dSBarry Smith   PetscFunctionBegin;
1723a40ed3dSBarry Smith   PetscFunctionReturn(0);
173ff14e315SSatish Balay }
174ff14e315SSatish Balay 
1754a2ae208SSatish Balay #undef __FUNCT__
1764a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
1778f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
1788965ea79SLois Curfman McInnes {
17939ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
1808965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
181d36fbae8SSatish Balay   int          ierr,nstash,reallocs;
1828965ea79SLois Curfman McInnes   InsertMode   addv;
1838965ea79SLois Curfman McInnes 
1843a40ed3dSBarry Smith   PetscFunctionBegin;
1858965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
186ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
1877056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
18829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
1898965ea79SLois Curfman McInnes   }
190e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
1918965ea79SLois Curfman McInnes 
1928798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
1938798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
194b0a32e0cSBarry Smith   PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs);
1953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1968965ea79SLois Curfman McInnes }
1978965ea79SLois Curfman McInnes 
1984a2ae208SSatish Balay #undef __FUNCT__
1994a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
2008f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2018965ea79SLois Curfman McInnes {
20239ddd567SLois Curfman McInnes   Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
2037ef1d9bdSSatish Balay   int          i,n,ierr,*row,*col,flg,j,rstart,ncols;
20487828ca2SBarry Smith   PetscScalar  *val;
205e0fa3b82SLois Curfman McInnes   InsertMode   addv=mat->insertmode;
2068965ea79SLois Curfman McInnes 
2073a40ed3dSBarry Smith   PetscFunctionBegin;
2088965ea79SLois Curfman McInnes   /*  wait on receives */
2097ef1d9bdSSatish Balay   while (1) {
2108798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2117ef1d9bdSSatish Balay     if (!flg) break;
2128965ea79SLois Curfman McInnes 
2137ef1d9bdSSatish Balay     for (i=0; i<n;) {
2147ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2157ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2167ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2177ef1d9bdSSatish Balay       else       ncols = n-i;
2187ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2197ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2207ef1d9bdSSatish Balay       i = j;
2218965ea79SLois Curfman McInnes     }
2227ef1d9bdSSatish Balay   }
2238798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2248965ea79SLois Curfman McInnes 
22539ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
22639ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2278965ea79SLois Curfman McInnes 
2286d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
22939ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2308965ea79SLois Curfman McInnes   }
2313a40ed3dSBarry Smith   PetscFunctionReturn(0);
2328965ea79SLois Curfman McInnes }
2338965ea79SLois Curfman McInnes 
2344a2ae208SSatish Balay #undef __FUNCT__
2354a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
2368f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A)
2378965ea79SLois Curfman McInnes {
2383a40ed3dSBarry Smith   int          ierr;
23939ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense*)A->data;
2403a40ed3dSBarry Smith 
2413a40ed3dSBarry Smith   PetscFunctionBegin;
2423a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2433a40ed3dSBarry Smith   PetscFunctionReturn(0);
2448965ea79SLois Curfman McInnes }
2458965ea79SLois Curfman McInnes 
2464a2ae208SSatish Balay #undef __FUNCT__
2474a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIDense"
2488f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs)
2494e220ebcSLois Curfman McInnes {
2503a40ed3dSBarry Smith   PetscFunctionBegin;
2514e220ebcSLois Curfman McInnes   *bs = 1;
2523a40ed3dSBarry Smith   PetscFunctionReturn(0);
2534e220ebcSLois Curfman McInnes }
2544e220ebcSLois Curfman McInnes 
2558965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2568965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2578965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2583501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2598965ea79SLois Curfman McInnes    routine.
2608965ea79SLois Curfman McInnes */
2614a2ae208SSatish Balay #undef __FUNCT__
2624a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
263268466fbSBarry Smith int MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag)
2648965ea79SLois Curfman McInnes {
26539ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
2668965ea79SLois Curfman McInnes   int            i,ierr,N,*rows,*owners = l->rowners,size = l->size;
267c1dc657dSBarry Smith   int            *nprocs,j,idx,nsends;
2688965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2698965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2708965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2718965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2728965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2738965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2748965ea79SLois Curfman McInnes   IS             istmp;
27535d8aa7fSBarry Smith   PetscTruth     found;
2768965ea79SLois Curfman McInnes 
2773a40ed3dSBarry Smith   PetscFunctionBegin;
278b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
2798965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
2808965ea79SLois Curfman McInnes 
2818965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
282b0a32e0cSBarry Smith   ierr  = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr);
283549d3d68SSatish Balay   ierr  = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
284b0a32e0cSBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/
2858965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
2868965ea79SLois Curfman McInnes     idx = rows[i];
28735d8aa7fSBarry Smith     found = PETSC_FALSE;
2888965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
2898965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
290c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
2918965ea79SLois Curfman McInnes       }
2928965ea79SLois Curfman McInnes     }
29329bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
2948965ea79SLois Curfman McInnes   }
295c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
2968965ea79SLois Curfman McInnes 
2978965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
298c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
2998965ea79SLois Curfman McInnes 
3008965ea79SLois Curfman McInnes   /* post receives:   */
301b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr);
302b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3038965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
304ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3058965ea79SLois Curfman McInnes   }
3068965ea79SLois Curfman McInnes 
3078965ea79SLois Curfman McInnes   /* do sends:
3088965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3098965ea79SLois Curfman McInnes          the ith processor
3108965ea79SLois Curfman McInnes   */
311b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr);
312b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
313b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr);
3148965ea79SLois Curfman McInnes   starts[0]  = 0;
315c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3168965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3178965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3188965ea79SLois Curfman McInnes   }
3198965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3208965ea79SLois Curfman McInnes 
3218965ea79SLois Curfman McInnes   starts[0] = 0;
322c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3238965ea79SLois Curfman McInnes   count = 0;
3248965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
325c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
326c1dc657dSBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3278965ea79SLois Curfman McInnes     }
3288965ea79SLois Curfman McInnes   }
329606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3308965ea79SLois Curfman McInnes 
3318965ea79SLois Curfman McInnes   base = owners[rank];
3328965ea79SLois Curfman McInnes 
3338965ea79SLois Curfman McInnes   /*  wait on receives */
334b0a32e0cSBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr);
3358965ea79SLois Curfman McInnes   source = lens + nrecvs;
3368965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3378965ea79SLois Curfman McInnes   while (count) {
338ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3398965ea79SLois Curfman McInnes     /* unpack receives into our local space */
340ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3418965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3428965ea79SLois Curfman McInnes     lens[imdex]  = n;
3438965ea79SLois Curfman McInnes     slen += n;
3448965ea79SLois Curfman McInnes     count--;
3458965ea79SLois Curfman McInnes   }
346606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3478965ea79SLois Curfman McInnes 
3488965ea79SLois Curfman McInnes   /* move the data into the send scatter */
349b0a32e0cSBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr);
3508965ea79SLois Curfman McInnes   count = 0;
3518965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
3528965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3538965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
3548965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3558965ea79SLois Curfman McInnes     }
3568965ea79SLois Curfman McInnes   }
357606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
358606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
359606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
360606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3618965ea79SLois Curfman McInnes 
3628965ea79SLois Curfman McInnes   /* actually zap the local rows */
363029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
364b0a32e0cSBarry Smith   PetscLogObjectParent(A,istmp);
365606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3668965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3678965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3688965ea79SLois Curfman McInnes 
3698965ea79SLois Curfman McInnes   /* wait on sends */
3708965ea79SLois Curfman McInnes   if (nsends) {
371b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
372ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
373606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
3748965ea79SLois Curfman McInnes   }
375606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
376606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes 
3783a40ed3dSBarry Smith   PetscFunctionReturn(0);
3798965ea79SLois Curfman McInnes }
3808965ea79SLois Curfman McInnes 
3814a2ae208SSatish Balay #undef __FUNCT__
3824a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
3838f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3848965ea79SLois Curfman McInnes {
38539ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
3868965ea79SLois Curfman McInnes   int          ierr;
387c456f294SBarry Smith 
3883a40ed3dSBarry Smith   PetscFunctionBegin;
38943a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39043a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39144cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
3923a40ed3dSBarry Smith   PetscFunctionReturn(0);
3938965ea79SLois Curfman McInnes }
3948965ea79SLois Curfman McInnes 
3954a2ae208SSatish Balay #undef __FUNCT__
3964a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
3978f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
3988965ea79SLois Curfman McInnes {
39939ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4008965ea79SLois Curfman McInnes   int          ierr;
401c456f294SBarry Smith 
4023a40ed3dSBarry Smith   PetscFunctionBegin;
40343a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40443a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40544cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4063a40ed3dSBarry Smith   PetscFunctionReturn(0);
4078965ea79SLois Curfman McInnes }
4088965ea79SLois Curfman McInnes 
4094a2ae208SSatish Balay #undef __FUNCT__
4104a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
4117c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
412096963f5SLois Curfman McInnes {
413096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
414096963f5SLois Curfman McInnes   int          ierr;
41587828ca2SBarry Smith   PetscScalar  zero = 0.0;
416096963f5SLois Curfman McInnes 
4173a40ed3dSBarry Smith   PetscFunctionBegin;
4183501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
4197c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
420537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
421537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4223a40ed3dSBarry Smith   PetscFunctionReturn(0);
423096963f5SLois Curfman McInnes }
424096963f5SLois Curfman McInnes 
4254a2ae208SSatish Balay #undef __FUNCT__
4264a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
4277c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
428096963f5SLois Curfman McInnes {
429096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
430096963f5SLois Curfman McInnes   int          ierr;
431096963f5SLois Curfman McInnes 
4323a40ed3dSBarry Smith   PetscFunctionBegin;
4333501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4347c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
435537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
436537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4373a40ed3dSBarry Smith   PetscFunctionReturn(0);
438096963f5SLois Curfman McInnes }
439096963f5SLois Curfman McInnes 
4404a2ae208SSatish Balay #undef __FUNCT__
4414a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
4428f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v)
4438965ea79SLois Curfman McInnes {
44439ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
445096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
446273d9f13SBarry Smith   int          ierr,len,i,n,m = A->m,radd;
44787828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
448ed3cc1f0SBarry Smith 
4493a40ed3dSBarry Smith   PetscFunctionBegin;
450273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
451*b1d4fb26SBarry Smith   ierr = VecGetArrayFast(v,&x);CHKERRQ(ierr);
452096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
453273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
454273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
4557ddc982cSLois Curfman McInnes   radd = a->rstart*m;
45644cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
457096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
458096963f5SLois Curfman McInnes   }
459*b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(v,&x);CHKERRQ(ierr);
4603a40ed3dSBarry Smith   PetscFunctionReturn(0);
4618965ea79SLois Curfman McInnes }
4628965ea79SLois Curfman McInnes 
4634a2ae208SSatish Balay #undef __FUNCT__
4644a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
465e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat)
4668965ea79SLois Curfman McInnes {
4673501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4688965ea79SLois Curfman McInnes   int          ierr;
469ed3cc1f0SBarry Smith 
4703a40ed3dSBarry Smith   PetscFunctionBegin;
47194d884c6SBarry Smith 
472aa482453SBarry Smith #if defined(PETSC_USE_LOG)
473b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N);
4748965ea79SLois Curfman McInnes #endif
4758798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
476606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
4773501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
4783501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4793501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
480622d7880SLois Curfman McInnes   if (mdn->factor) {
481606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
482606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
483606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
484606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
485622d7880SLois Curfman McInnes   }
486606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
4873a40ed3dSBarry Smith   PetscFunctionReturn(0);
4888965ea79SLois Curfman McInnes }
48939ddd567SLois Curfman McInnes 
4904a2ae208SSatish Balay #undef __FUNCT__
4914a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
492b0a32e0cSBarry Smith static int MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
4938965ea79SLois Curfman McInnes {
49439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4958965ea79SLois Curfman McInnes   int          ierr;
4967056b6fcSBarry Smith 
4973a40ed3dSBarry Smith   PetscFunctionBegin;
49839ddd567SLois Curfman McInnes   if (mdn->size == 1) {
49939ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5008965ea79SLois Curfman McInnes   }
50129bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5023a40ed3dSBarry Smith   PetscFunctionReturn(0);
5038965ea79SLois Curfman McInnes }
5048965ea79SLois Curfman McInnes 
5054a2ae208SSatish Balay #undef __FUNCT__
5064a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
507b0a32e0cSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5088965ea79SLois Curfman McInnes {
50939ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
510fb9695e5SSatish Balay   int               ierr,size = mdn->size,rank = mdn->rank;
511b0a32e0cSBarry Smith   PetscViewerType   vtype;
512f1af5d2fSBarry Smith   PetscTruth        isascii,isdraw;
513b0a32e0cSBarry Smith   PetscViewer       sviewer;
514f3ef73ceSBarry Smith   PetscViewerFormat format;
5158965ea79SLois Curfman McInnes 
5163a40ed3dSBarry Smith   PetscFunctionBegin;
517b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
518fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
519f1af5d2fSBarry Smith   if (isascii) {
520b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
521b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
522456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
5234e220ebcSLois Curfman McInnes       MatInfo info;
524888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
525b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m,
5266831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
527b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5283501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5293a40ed3dSBarry Smith       PetscFunctionReturn(0);
530fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5313a40ed3dSBarry Smith       PetscFunctionReturn(0);
5328965ea79SLois Curfman McInnes     }
533f1af5d2fSBarry Smith   } else if (isdraw) {
534b0a32e0cSBarry Smith     PetscDraw       draw;
535f1af5d2fSBarry Smith     PetscTruth isnull;
536f1af5d2fSBarry Smith 
537b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
538b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
539f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
540f1af5d2fSBarry Smith   }
54177ed5343SBarry Smith 
5428965ea79SLois Curfman McInnes   if (size == 1) {
54339ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5443a40ed3dSBarry Smith   } else {
5458965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5468965ea79SLois Curfman McInnes     Mat          A;
547273d9f13SBarry Smith     int          M = mat->M,N = mat->N,m,row,i,nz,*cols;
54887828ca2SBarry Smith     PetscScalar  *vals;
5498965ea79SLois Curfman McInnes 
5508965ea79SLois Curfman McInnes     if (!rank) {
551f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5523a40ed3dSBarry Smith     } else {
553f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5548965ea79SLois Curfman McInnes     }
555b0a32e0cSBarry Smith     PetscLogObjectParent(mat,A);
5568965ea79SLois Curfman McInnes 
55739ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
55839ddd567SLois Curfman McInnes        but it's quick for now */
559273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
5608965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
56139ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
56239ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
56339ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
56439ddd567SLois Curfman McInnes       row++;
5658965ea79SLois Curfman McInnes     }
5668965ea79SLois Curfman McInnes 
5676d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5686d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
569b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
570b9b97703SBarry Smith     if (!rank) {
5716831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5728965ea79SLois Curfman McInnes     }
573b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
574b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5758965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
5768965ea79SLois Curfman McInnes   }
5773a40ed3dSBarry Smith   PetscFunctionReturn(0);
5788965ea79SLois Curfman McInnes }
5798965ea79SLois Curfman McInnes 
5804a2ae208SSatish Balay #undef __FUNCT__
5814a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
582b0a32e0cSBarry Smith int MatView_MPIDense(Mat mat,PetscViewer viewer)
5838965ea79SLois Curfman McInnes {
58439ddd567SLois Curfman McInnes   int        ierr;
585f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
5868965ea79SLois Curfman McInnes 
587433994e6SBarry Smith   PetscFunctionBegin;
5880f5bd95cSBarry Smith 
589b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
590fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
591b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
592fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
5930f5bd95cSBarry Smith 
594f1af5d2fSBarry Smith   if (isascii || issocket || isdraw) {
595f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
5960f5bd95cSBarry Smith   } else if (isbinary) {
5973a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
5985cd90555SBarry Smith   } else {
59929bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6008965ea79SLois Curfman McInnes   }
6013a40ed3dSBarry Smith   PetscFunctionReturn(0);
6028965ea79SLois Curfman McInnes }
6038965ea79SLois Curfman McInnes 
6044a2ae208SSatish Balay #undef __FUNCT__
6054a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
6068f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6078965ea79SLois Curfman McInnes {
6083501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6093501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6104e220ebcSLois Curfman McInnes   int          ierr;
611329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6128965ea79SLois Curfman McInnes 
6133a40ed3dSBarry Smith   PetscFunctionBegin;
614273d9f13SBarry Smith   info->rows_global    = (double)A->M;
615273d9f13SBarry Smith   info->columns_global = (double)A->N;
616273d9f13SBarry Smith   info->rows_local     = (double)A->m;
617273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6184e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6194e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6204e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6214e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6228965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6234e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6244e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6254e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6264e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6274e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6288965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
629d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
6304e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6314e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6324e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6334e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6344e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6358965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
636d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
6374e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6384e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6394e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6404e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6414e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6428965ea79SLois Curfman McInnes   }
6434e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6444e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6454e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6463a40ed3dSBarry Smith   PetscFunctionReturn(0);
6478965ea79SLois Curfman McInnes }
6488965ea79SLois Curfman McInnes 
6494a2ae208SSatish Balay #undef __FUNCT__
6504a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
6518f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6528965ea79SLois Curfman McInnes {
65339ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
654273d9f13SBarry Smith   int          ierr;
6558965ea79SLois Curfman McInnes 
6563a40ed3dSBarry Smith   PetscFunctionBegin;
65712c028f9SKris Buschelman   switch (op) {
65812c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
65912c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
66012c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
66112c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
66212c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
66312c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
664273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
66512c028f9SKris Buschelman     break;
66612c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
667273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
668273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
66912c028f9SKris Buschelman     break;
67012c028f9SKris Buschelman   case MAT_ROWS_SORTED:
67112c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
67212c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
67312c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
674b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
67512c028f9SKris Buschelman     break;
67612c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
677273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
678273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
67912c028f9SKris Buschelman     break;
68012c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
681273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
68212c028f9SKris Buschelman     break;
68312c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
68429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
68577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
68677e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
6879a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
6889a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
6899a4540c5SBarry Smith   case MAT_HERMITIAN:
6909a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
6919a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
6929a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
69377e54ba9SKris Buschelman     break;
69412c028f9SKris Buschelman   default:
69529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
6963a40ed3dSBarry Smith   }
6973a40ed3dSBarry Smith   PetscFunctionReturn(0);
6988965ea79SLois Curfman McInnes }
6998965ea79SLois Curfman McInnes 
7004a2ae208SSatish Balay #undef __FUNCT__
7014a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIDense"
70287828ca2SBarry Smith int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v)
7038965ea79SLois Curfman McInnes {
7043501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7053a40ed3dSBarry Smith   int          lrow,rstart = mat->rstart,rend = mat->rend,ierr;
7068965ea79SLois Curfman McInnes 
7073a40ed3dSBarry Smith   PetscFunctionBegin;
70829bbc08cSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
7098965ea79SLois Curfman McInnes   lrow = row - rstart;
7103a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7113a40ed3dSBarry Smith   PetscFunctionReturn(0);
7128965ea79SLois Curfman McInnes }
7138965ea79SLois Curfman McInnes 
7144a2ae208SSatish Balay #undef __FUNCT__
7154a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIDense"
71687828ca2SBarry Smith int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v)
7178965ea79SLois Curfman McInnes {
718606d414cSSatish Balay   int ierr;
719606d414cSSatish Balay 
7203a40ed3dSBarry Smith   PetscFunctionBegin;
721606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
722606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7233a40ed3dSBarry Smith   PetscFunctionReturn(0);
7248965ea79SLois Curfman McInnes }
7258965ea79SLois Curfman McInnes 
7264a2ae208SSatish Balay #undef __FUNCT__
7274a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
7285b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7295b2fa520SLois Curfman McInnes {
7305b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7315b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
73287828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
733273d9f13SBarry Smith   int          ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7345b2fa520SLois Curfman McInnes 
7355b2fa520SLois Curfman McInnes   PetscFunctionBegin;
73672d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7375b2fa520SLois Curfman McInnes   if (ll) {
73872d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
73929bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
740*b1d4fb26SBarry Smith     ierr = VecGetArrayFast(ll,&l);CHKERRQ(ierr);
7415b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7425b2fa520SLois Curfman McInnes       x = l[i];
7435b2fa520SLois Curfman McInnes       v = mat->v + i;
7445b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7455b2fa520SLois Curfman McInnes     }
746*b1d4fb26SBarry Smith     ierr = VecRestoreArrayFast(ll,&l);CHKERRQ(ierr);
747b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7485b2fa520SLois Curfman McInnes   }
7495b2fa520SLois Curfman McInnes   if (rr) {
75072d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
75129bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7525b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7535b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
754*b1d4fb26SBarry Smith     ierr = VecGetArrayFast(mdn->lvec,&r);CHKERRQ(ierr);
7555b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7565b2fa520SLois Curfman McInnes       x = r[i];
7575b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7585b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7595b2fa520SLois Curfman McInnes     }
760*b1d4fb26SBarry Smith     ierr = VecRestoreArrayFast(mdn->lvec,&r);CHKERRQ(ierr);
761b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7625b2fa520SLois Curfman McInnes   }
7635b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7645b2fa520SLois Curfman McInnes }
7655b2fa520SLois Curfman McInnes 
7664a2ae208SSatish Balay #undef __FUNCT__
7674a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
768064f8208SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
769096963f5SLois Curfman McInnes {
7703501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7713501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7723501a2bdSLois Curfman McInnes   int          ierr,i,j;
773329f5518SBarry Smith   PetscReal    sum = 0.0;
77487828ca2SBarry Smith   PetscScalar  *v = mat->v;
7753501a2bdSLois Curfman McInnes 
7763a40ed3dSBarry Smith   PetscFunctionBegin;
7773501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
778064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
7793501a2bdSLois Curfman McInnes   } else {
7803501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
781273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
782aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
783329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
7843501a2bdSLois Curfman McInnes #else
7853501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7863501a2bdSLois Curfman McInnes #endif
7873501a2bdSLois Curfman McInnes       }
788064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
789064f8208SBarry Smith       *nrm = sqrt(*nrm);
790b0a32e0cSBarry Smith       PetscLogFlops(2*mdn->A->n*mdn->A->m);
7913a40ed3dSBarry Smith     } else if (type == NORM_1) {
792329f5518SBarry Smith       PetscReal *tmp,*tmp2;
793b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
794273d9f13SBarry Smith       tmp2 = tmp + A->N;
795273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
796064f8208SBarry Smith       *nrm = 0.0;
7973501a2bdSLois Curfman McInnes       v = mat->v;
798273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
799273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
80067e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8013501a2bdSLois Curfman McInnes         }
8023501a2bdSLois Curfman McInnes       }
803d7d1e502SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
804273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
805064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8063501a2bdSLois Curfman McInnes       }
807606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
808b0a32e0cSBarry Smith       PetscLogFlops(A->n*A->m);
8093a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
810329f5518SBarry Smith       PetscReal ntemp;
8113501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
812064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
8133a40ed3dSBarry Smith     } else {
81429bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8153501a2bdSLois Curfman McInnes     }
8163501a2bdSLois Curfman McInnes   }
8173a40ed3dSBarry Smith   PetscFunctionReturn(0);
8183501a2bdSLois Curfman McInnes }
8193501a2bdSLois Curfman McInnes 
8204a2ae208SSatish Balay #undef __FUNCT__
8214a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
8228f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8233501a2bdSLois Curfman McInnes {
8243501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8253501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8263501a2bdSLois Curfman McInnes   Mat          B;
827273d9f13SBarry Smith   int          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8283501a2bdSLois Curfman McInnes   int          j,i,ierr;
82987828ca2SBarry Smith   PetscScalar  *v;
8303501a2bdSLois Curfman McInnes 
8313a40ed3dSBarry Smith   PetscFunctionBegin;
8327c922b88SBarry Smith   if (!matout && M != N) {
83329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8347056b6fcSBarry Smith   }
8357056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8363501a2bdSLois Curfman McInnes 
837273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
838b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr);
8393501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8403501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8413501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8423501a2bdSLois Curfman McInnes     v   += m;
8433501a2bdSLois Curfman McInnes   }
844606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8456d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8466d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8477c922b88SBarry Smith   if (matout) {
8483501a2bdSLois Curfman McInnes     *matout = B;
8493501a2bdSLois Curfman McInnes   } else {
850273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8513501a2bdSLois Curfman McInnes   }
8523a40ed3dSBarry Smith   PetscFunctionReturn(0);
853096963f5SLois Curfman McInnes }
854096963f5SLois Curfman McInnes 
855d9eff348SSatish Balay #include "petscblaslapack.h"
8564a2ae208SSatish Balay #undef __FUNCT__
8574a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
858268466fbSBarry Smith int MatScale_MPIDense(const PetscScalar *alpha,Mat inA)
85944cd7ae7SLois Curfman McInnes {
86044cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
86144cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
86244cd7ae7SLois Curfman McInnes   int          one = 1,nz;
86344cd7ae7SLois Curfman McInnes 
8643a40ed3dSBarry Smith   PetscFunctionBegin;
865273d9f13SBarry Smith   nz = inA->m*inA->N;
866268466fbSBarry Smith   BLscal_(&nz,(PetscScalar*)alpha,a->v,&one);
867b0a32e0cSBarry Smith   PetscLogFlops(nz);
8683a40ed3dSBarry Smith   PetscFunctionReturn(0);
86944cd7ae7SLois Curfman McInnes }
87044cd7ae7SLois Curfman McInnes 
8715609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
8728965ea79SLois Curfman McInnes 
8734a2ae208SSatish Balay #undef __FUNCT__
8744a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
875273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A)
876273d9f13SBarry Smith {
877273d9f13SBarry Smith   int        ierr;
878273d9f13SBarry Smith 
879273d9f13SBarry Smith   PetscFunctionBegin;
880273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
881273d9f13SBarry Smith   PetscFunctionReturn(0);
882273d9f13SBarry Smith }
883273d9f13SBarry Smith 
8848965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
88509dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
88609dc0095SBarry Smith        MatGetRow_MPIDense,
88709dc0095SBarry Smith        MatRestoreRow_MPIDense,
88809dc0095SBarry Smith        MatMult_MPIDense,
88997304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
8907c922b88SBarry Smith        MatMultTranspose_MPIDense,
8917c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
8928965ea79SLois Curfman McInnes        0,
89309dc0095SBarry Smith        0,
89409dc0095SBarry Smith        0,
89597304618SKris Buschelman /*10*/ 0,
89609dc0095SBarry Smith        0,
89709dc0095SBarry Smith        0,
89809dc0095SBarry Smith        0,
89909dc0095SBarry Smith        MatTranspose_MPIDense,
90097304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
90197304618SKris Buschelman        0,
90209dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9035b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
90409dc0095SBarry Smith        MatNorm_MPIDense,
90597304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
90609dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
90709dc0095SBarry Smith        0,
90809dc0095SBarry Smith        MatSetOption_MPIDense,
90909dc0095SBarry Smith        MatZeroEntries_MPIDense,
91097304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
91109dc0095SBarry Smith        0,
91209dc0095SBarry Smith        0,
91309dc0095SBarry Smith        0,
91409dc0095SBarry Smith        0,
91597304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
916273d9f13SBarry Smith        0,
91709dc0095SBarry Smith        0,
91809dc0095SBarry Smith        MatGetArray_MPIDense,
91909dc0095SBarry Smith        MatRestoreArray_MPIDense,
92097304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
92109dc0095SBarry Smith        0,
92209dc0095SBarry Smith        0,
92309dc0095SBarry Smith        0,
92409dc0095SBarry Smith        0,
92597304618SKris Buschelman /*40*/ 0,
9262ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
92709dc0095SBarry Smith        0,
92809dc0095SBarry Smith        MatGetValues_MPIDense,
92909dc0095SBarry Smith        0,
93097304618SKris Buschelman /*45*/ 0,
93109dc0095SBarry Smith        MatScale_MPIDense,
93209dc0095SBarry Smith        0,
93309dc0095SBarry Smith        0,
93409dc0095SBarry Smith        0,
93597304618SKris Buschelman /*50*/ MatGetBlockSize_MPIDense,
93609dc0095SBarry Smith        0,
93709dc0095SBarry Smith        0,
93809dc0095SBarry Smith        0,
93909dc0095SBarry Smith        0,
94097304618SKris Buschelman /*55*/ 0,
94109dc0095SBarry Smith        0,
94209dc0095SBarry Smith        0,
94309dc0095SBarry Smith        0,
94409dc0095SBarry Smith        0,
94597304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
946b9b97703SBarry Smith        MatDestroy_MPIDense,
947b9b97703SBarry Smith        MatView_MPIDense,
94897304618SKris Buschelman        MatGetPetscMaps_Petsc,
94997304618SKris Buschelman        0,
95097304618SKris Buschelman /*65*/ 0,
95197304618SKris Buschelman        0,
95297304618SKris Buschelman        0,
95397304618SKris Buschelman        0,
95497304618SKris Buschelman        0,
95597304618SKris Buschelman /*70*/ 0,
95697304618SKris Buschelman        0,
95797304618SKris Buschelman        0,
95897304618SKris Buschelman        0,
95997304618SKris Buschelman        0,
96097304618SKris Buschelman /*75*/ 0,
96197304618SKris Buschelman        0,
96297304618SKris Buschelman        0,
96397304618SKris Buschelman        0,
96497304618SKris Buschelman        0,
96597304618SKris Buschelman /*80*/ 0,
96697304618SKris Buschelman        0,
96797304618SKris Buschelman        0,
96897304618SKris Buschelman        0,
96997304618SKris Buschelman /*85*/ MatLoad_MPIDense};
9708965ea79SLois Curfman McInnes 
971273d9f13SBarry Smith EXTERN_C_BEGIN
9724a2ae208SSatish Balay #undef __FUNCT__
973a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
974a23d5eceSKris Buschelman int MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
975a23d5eceSKris Buschelman {
976a23d5eceSKris Buschelman   Mat_MPIDense *a;
977a23d5eceSKris Buschelman   int          ierr;
978a23d5eceSKris Buschelman 
979a23d5eceSKris Buschelman   PetscFunctionBegin;
980a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
981a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
982a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
983a23d5eceSKris Buschelman 
984a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
985a23d5eceSKris Buschelman   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mat->m,mat->N,data,&a->A);CHKERRQ(ierr);
986a23d5eceSKris Buschelman   PetscLogObjectParent(mat,a->A);
987a23d5eceSKris Buschelman   PetscFunctionReturn(0);
988a23d5eceSKris Buschelman }
989a23d5eceSKris Buschelman EXTERN_C_END
990a23d5eceSKris Buschelman 
9910bad9183SKris Buschelman /*MC
992fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
9930bad9183SKris Buschelman 
9940bad9183SKris Buschelman    Options Database Keys:
9950bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
9960bad9183SKris Buschelman 
9970bad9183SKris Buschelman   Level: beginner
9980bad9183SKris Buschelman 
9990bad9183SKris Buschelman .seealso: MatCreateMPIDense
10000bad9183SKris Buschelman M*/
10010bad9183SKris Buschelman 
1002a23d5eceSKris Buschelman EXTERN_C_BEGIN
1003a23d5eceSKris Buschelman #undef __FUNCT__
10044a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1005273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat)
1006273d9f13SBarry Smith {
1007273d9f13SBarry Smith   Mat_MPIDense *a;
1008273d9f13SBarry Smith   int          ierr,i;
1009273d9f13SBarry Smith 
1010273d9f13SBarry Smith   PetscFunctionBegin;
1011b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1012b0a32e0cSBarry Smith   mat->data         = (void*)a;
1013273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1014273d9f13SBarry Smith   mat->factor       = 0;
1015273d9f13SBarry Smith   mat->mapping      = 0;
1016273d9f13SBarry Smith 
1017273d9f13SBarry Smith   a->factor       = 0;
1018273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1019273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
1020273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
1021273d9f13SBarry Smith 
1022273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
1023273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
1024273d9f13SBarry Smith   a->nvec = mat->n;
1025273d9f13SBarry Smith 
1026273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1027273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
10288a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
10298a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
1030273d9f13SBarry Smith 
1031273d9f13SBarry Smith   /* build local table of row and column ownerships */
1032b0a32e0cSBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1033273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
1034b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1035273d9f13SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1036273d9f13SBarry Smith   a->rowners[0] = 0;
1037273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1038273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
1039273d9f13SBarry Smith   }
1040273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
1041273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
1042273d9f13SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1043273d9f13SBarry Smith   a->cowners[0] = 0;
1044273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1045273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1046273d9f13SBarry Smith   }
1047273d9f13SBarry Smith 
1048273d9f13SBarry Smith   /* build cache for off array entries formed */
1049273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1050273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1051273d9f13SBarry Smith 
1052273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1053273d9f13SBarry Smith   a->lvec        = 0;
1054273d9f13SBarry Smith   a->Mvctx       = 0;
1055273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1056273d9f13SBarry Smith 
1057273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1058273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1059273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1060a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1061a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1062a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1063273d9f13SBarry Smith   PetscFunctionReturn(0);
1064273d9f13SBarry Smith }
1065273d9f13SBarry Smith EXTERN_C_END
1066273d9f13SBarry Smith 
1067209238afSKris Buschelman /*MC
1068002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1069209238afSKris Buschelman 
1070209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1071209238afSKris Buschelman    and MATMPIDENSE otherwise.
1072209238afSKris Buschelman 
1073209238afSKris Buschelman    Options Database Keys:
1074209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1075209238afSKris Buschelman 
1076209238afSKris Buschelman   Level: beginner
1077209238afSKris Buschelman 
1078209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1079209238afSKris Buschelman M*/
1080209238afSKris Buschelman 
1081209238afSKris Buschelman EXTERN_C_BEGIN
1082209238afSKris Buschelman #undef __FUNCT__
1083209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1084209238afSKris Buschelman int MatCreate_Dense(Mat A) {
1085209238afSKris Buschelman   int ierr,size;
1086209238afSKris Buschelman 
1087209238afSKris Buschelman   PetscFunctionBegin;
1088209238afSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr);
1089209238afSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1090209238afSKris Buschelman   if (size == 1) {
1091209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1092209238afSKris Buschelman   } else {
1093209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1094209238afSKris Buschelman   }
1095209238afSKris Buschelman   PetscFunctionReturn(0);
1096209238afSKris Buschelman }
1097209238afSKris Buschelman EXTERN_C_END
1098209238afSKris Buschelman 
10994a2ae208SSatish Balay #undef __FUNCT__
11004a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1101273d9f13SBarry Smith /*@C
1102273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1103273d9f13SBarry Smith 
1104273d9f13SBarry Smith    Not collective
1105273d9f13SBarry Smith 
1106273d9f13SBarry Smith    Input Parameters:
1107273d9f13SBarry Smith .  A - the matrix
1108273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1109273d9f13SBarry Smith    to control all matrix memory allocation.
1110273d9f13SBarry Smith 
1111273d9f13SBarry Smith    Notes:
1112273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1113273d9f13SBarry Smith    storage by columns.
1114273d9f13SBarry Smith 
1115273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1116273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1117273d9f13SBarry Smith    set data=PETSC_NULL.
1118273d9f13SBarry Smith 
1119273d9f13SBarry Smith    Level: intermediate
1120273d9f13SBarry Smith 
1121273d9f13SBarry Smith .keywords: matrix,dense, parallel
1122273d9f13SBarry Smith 
1123273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1124273d9f13SBarry Smith @*/
112587828ca2SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1126273d9f13SBarry Smith {
1127a23d5eceSKris Buschelman   int ierr,(*f)(Mat,PetscScalar *);
1128273d9f13SBarry Smith 
1129273d9f13SBarry Smith   PetscFunctionBegin;
1130565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1131a23d5eceSKris Buschelman   if (f) {
1132a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1133a23d5eceSKris Buschelman   }
1134273d9f13SBarry Smith   PetscFunctionReturn(0);
1135273d9f13SBarry Smith }
1136273d9f13SBarry Smith 
11374a2ae208SSatish Balay #undef __FUNCT__
11384a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
11398965ea79SLois Curfman McInnes /*@C
114039ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
11418965ea79SLois Curfman McInnes 
1142db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1143db81eaa0SLois Curfman McInnes 
11448965ea79SLois Curfman McInnes    Input Parameters:
1145db81eaa0SLois Curfman McInnes +  comm - MPI communicator
11468965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1147db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
11488965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1149db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
11507f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1151dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
11528965ea79SLois Curfman McInnes 
11538965ea79SLois Curfman McInnes    Output Parameter:
1154477f1c0bSLois Curfman McInnes .  A - the matrix
11558965ea79SLois Curfman McInnes 
1156b259b22eSLois Curfman McInnes    Notes:
115739ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
115839ddd567SLois Curfman McInnes    storage by columns.
11598965ea79SLois Curfman McInnes 
116018f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
116118f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
11627f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
116318f449edSLois Curfman McInnes 
11648965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
11658965ea79SLois Curfman McInnes    (possibly both).
11668965ea79SLois Curfman McInnes 
1167027ccd11SLois Curfman McInnes    Level: intermediate
1168027ccd11SLois Curfman McInnes 
116939ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
11708965ea79SLois Curfman McInnes 
117139ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
11728965ea79SLois Curfman McInnes @*/
117387828ca2SBarry Smith int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A)
11748965ea79SLois Curfman McInnes {
1175273d9f13SBarry Smith   int ierr,size;
11768965ea79SLois Curfman McInnes 
11773a40ed3dSBarry Smith   PetscFunctionBegin;
1178273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1179273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1180273d9f13SBarry Smith   if (size > 1) {
1181273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1182273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1183273d9f13SBarry Smith   } else {
1184273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1185273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
11868c469469SLois Curfman McInnes   }
11873a40ed3dSBarry Smith   PetscFunctionReturn(0);
11888965ea79SLois Curfman McInnes }
11898965ea79SLois Curfman McInnes 
11904a2ae208SSatish Balay #undef __FUNCT__
11914a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
11925609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11938965ea79SLois Curfman McInnes {
11948965ea79SLois Curfman McInnes   Mat          mat;
11953501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
119639ddd567SLois Curfman McInnes   int          ierr;
11978965ea79SLois Curfman McInnes 
11983a40ed3dSBarry Smith   PetscFunctionBegin;
11998965ea79SLois Curfman McInnes   *newmat       = 0;
1200273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1201273d9f13SBarry Smith   ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr);
1202b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1203b0a32e0cSBarry Smith   mat->data         = (void*)a;
1204549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
12053501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1206c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1207273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
12088965ea79SLois Curfman McInnes 
12098965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
12108965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
12118965ea79SLois Curfman McInnes   a->size         = oldmat->size;
12128965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1213e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1214b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
12153782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1216b0a32e0cSBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1217b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1218549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
12198798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
12208965ea79SLois Curfman McInnes 
1221329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
12225609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
1223b0a32e0cSBarry Smith   PetscLogObjectParent(mat,a->A);
12248965ea79SLois Curfman McInnes   *newmat = mat;
12253a40ed3dSBarry Smith   PetscFunctionReturn(0);
12268965ea79SLois Curfman McInnes }
12278965ea79SLois Curfman McInnes 
1228e090d566SSatish Balay #include "petscsys.h"
12298965ea79SLois Curfman McInnes 
12304a2ae208SSatish Balay #undef __FUNCT__
12314a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
123290ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat)
123390ace30eSBarry Smith {
123440011551SBarry Smith   int          *rowners,i,size,rank,m,ierr,nz,j;
123587828ca2SBarry Smith   PetscScalar  *array,*vals,*vals_ptr;
123690ace30eSBarry Smith   MPI_Status   status;
123790ace30eSBarry Smith 
12383a40ed3dSBarry Smith   PetscFunctionBegin;
1239d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1240d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
124190ace30eSBarry Smith 
124290ace30eSBarry Smith   /* determine ownership of all rows */
124390ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
1244b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1245ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
124690ace30eSBarry Smith   rowners[0] = 0;
124790ace30eSBarry Smith   for (i=2; i<=size; i++) {
124890ace30eSBarry Smith     rowners[i] += rowners[i-1];
124990ace30eSBarry Smith   }
125090ace30eSBarry Smith 
125190ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
125290ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
125390ace30eSBarry Smith 
125490ace30eSBarry Smith   if (!rank) {
125587828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
125690ace30eSBarry Smith 
125790ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
12580752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
125990ace30eSBarry Smith 
126090ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
126190ace30eSBarry Smith     vals_ptr = vals;
126290ace30eSBarry Smith     for (i=0; i<m; i++) {
126390ace30eSBarry Smith       for (j=0; j<N; j++) {
126490ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
126590ace30eSBarry Smith       }
126690ace30eSBarry Smith     }
126790ace30eSBarry Smith 
126890ace30eSBarry Smith     /* read in other processors and ship out */
126990ace30eSBarry Smith     for (i=1; i<size; i++) {
127090ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12710752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1272ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
127390ace30eSBarry Smith     }
12743a40ed3dSBarry Smith   } else {
127590ace30eSBarry Smith     /* receive numeric values */
127687828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
127790ace30eSBarry Smith 
127890ace30eSBarry Smith     /* receive message of values*/
1279ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
128090ace30eSBarry Smith 
128190ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
128290ace30eSBarry Smith     vals_ptr = vals;
128390ace30eSBarry Smith     for (i=0; i<m; i++) {
128490ace30eSBarry Smith       for (j=0; j<N; j++) {
128590ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
128690ace30eSBarry Smith       }
128790ace30eSBarry Smith     }
128890ace30eSBarry Smith   }
1289606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1290606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
12916d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12926d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12933a40ed3dSBarry Smith   PetscFunctionReturn(0);
129490ace30eSBarry Smith }
129590ace30eSBarry Smith 
12964a2ae208SSatish Balay #undef __FUNCT__
12974a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
12988e9aea5cSBarry Smith int MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat)
12998965ea79SLois Curfman McInnes {
13008965ea79SLois Curfman McInnes   Mat          A;
130187828ca2SBarry Smith   PetscScalar  *vals,*svals;
130219bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
13038965ea79SLois Curfman McInnes   MPI_Status   status;
13048965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
13058965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols;
130619bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
13073a40ed3dSBarry Smith   int          i,nz,ierr,j,rstart,rend,fd;
13088965ea79SLois Curfman McInnes 
13093a40ed3dSBarry Smith   PetscFunctionBegin;
1310d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1311d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13128965ea79SLois Curfman McInnes   if (!rank) {
1313b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
13140752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1315552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
13168965ea79SLois Curfman McInnes   }
13178965ea79SLois Curfman McInnes 
1318ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
131990ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
132090ace30eSBarry Smith 
132190ace30eSBarry Smith   /*
132290ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
132390ace30eSBarry Smith   */
132490ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
13253a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
13263a40ed3dSBarry Smith     PetscFunctionReturn(0);
132790ace30eSBarry Smith   }
132890ace30eSBarry Smith 
13298965ea79SLois Curfman McInnes   /* determine ownership of all rows */
13308965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
1331b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1332ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
13338965ea79SLois Curfman McInnes   rowners[0] = 0;
13348965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
13358965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
13368965ea79SLois Curfman McInnes   }
13378965ea79SLois Curfman McInnes   rstart = rowners[rank];
13388965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
13398965ea79SLois Curfman McInnes 
13408965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
1341b0a32e0cSBarry Smith   ierr    = PetscMalloc(2*(rend-rstart)*sizeof(int),&ourlens);CHKERRQ(ierr);
13428965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
13438965ea79SLois Curfman McInnes   if (!rank) {
1344b0a32e0cSBarry Smith     ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr);
13450752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1346b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr);
13478965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1348ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1349606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1350ca161407SBarry Smith   } else {
1351ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
13528965ea79SLois Curfman McInnes   }
13538965ea79SLois Curfman McInnes 
13548965ea79SLois Curfman McInnes   if (!rank) {
13558965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
1356b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr);
1357549d3d68SSatish Balay     ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
13588965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13598965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
13608965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
13618965ea79SLois Curfman McInnes       }
13628965ea79SLois Curfman McInnes     }
1363606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
13648965ea79SLois Curfman McInnes 
13658965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
13668965ea79SLois Curfman McInnes     maxnz = 0;
13678965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13680452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
13698965ea79SLois Curfman McInnes     }
1370b0a32e0cSBarry Smith     ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr);
13718965ea79SLois Curfman McInnes 
13728965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13738965ea79SLois Curfman McInnes     nz = procsnz[0];
1374b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13750752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13768965ea79SLois Curfman McInnes 
13778965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13788965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13798965ea79SLois Curfman McInnes       nz   = procsnz[i];
13800752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1381ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13828965ea79SLois Curfman McInnes     }
1383606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13843a40ed3dSBarry Smith   } else {
13858965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13868965ea79SLois Curfman McInnes     nz = 0;
13878965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13888965ea79SLois Curfman McInnes       nz += ourlens[i];
13898965ea79SLois Curfman McInnes     }
1390b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13918965ea79SLois Curfman McInnes 
13928965ea79SLois Curfman McInnes     /* receive message of column indices*/
1393ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1394ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
139529bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13968965ea79SLois Curfman McInnes   }
13978965ea79SLois Curfman McInnes 
13988965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1399549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
14008965ea79SLois Curfman McInnes   jj = 0;
14018965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14028965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
14038965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
14048965ea79SLois Curfman McInnes       jj++;
14058965ea79SLois Curfman McInnes     }
14068965ea79SLois Curfman McInnes   }
14078965ea79SLois Curfman McInnes 
14088965ea79SLois Curfman McInnes   /* create our matrix */
14098965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14108965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
14118965ea79SLois Curfman McInnes   }
1412b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
14138965ea79SLois Curfman McInnes   A = *newmat;
14148965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14158965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
14168965ea79SLois Curfman McInnes   }
14178965ea79SLois Curfman McInnes 
14188965ea79SLois Curfman McInnes   if (!rank) {
141987828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14208965ea79SLois Curfman McInnes 
14218965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
14228965ea79SLois Curfman McInnes     nz = procsnz[0];
14230752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
14248965ea79SLois Curfman McInnes 
14258965ea79SLois Curfman McInnes     /* insert into matrix */
14268965ea79SLois Curfman McInnes     jj      = rstart;
14278965ea79SLois Curfman McInnes     smycols = mycols;
14288965ea79SLois Curfman McInnes     svals   = vals;
14298965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14308965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14318965ea79SLois Curfman McInnes       smycols += ourlens[i];
14328965ea79SLois Curfman McInnes       svals   += ourlens[i];
14338965ea79SLois Curfman McInnes       jj++;
14348965ea79SLois Curfman McInnes     }
14358965ea79SLois Curfman McInnes 
14368965ea79SLois Curfman McInnes     /* read in other processors and ship out */
14378965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14388965ea79SLois Curfman McInnes       nz   = procsnz[i];
14390752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1440ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
14418965ea79SLois Curfman McInnes     }
1442606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
14433a40ed3dSBarry Smith   } else {
14448965ea79SLois Curfman McInnes     /* receive numeric values */
144587828ca2SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14468965ea79SLois Curfman McInnes 
14478965ea79SLois Curfman McInnes     /* receive message of values*/
1448ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1449ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
145029bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14518965ea79SLois Curfman McInnes 
14528965ea79SLois Curfman McInnes     /* insert into matrix */
14538965ea79SLois Curfman McInnes     jj      = rstart;
14548965ea79SLois Curfman McInnes     smycols = mycols;
14558965ea79SLois Curfman McInnes     svals   = vals;
14568965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14578965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14588965ea79SLois Curfman McInnes       smycols += ourlens[i];
14598965ea79SLois Curfman McInnes       svals   += ourlens[i];
14608965ea79SLois Curfman McInnes       jj++;
14618965ea79SLois Curfman McInnes     }
14628965ea79SLois Curfman McInnes   }
1463606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1464606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1465606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1466606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
14678965ea79SLois Curfman McInnes 
14686d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14696d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14703a40ed3dSBarry Smith   PetscFunctionReturn(0);
14718965ea79SLois Curfman McInnes }
147290ace30eSBarry Smith 
147390ace30eSBarry Smith 
147490ace30eSBarry Smith 
147590ace30eSBarry Smith 
147690ace30eSBarry Smith 
1477