xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision f3ef73ce49bf99b5a752264ab8d9443243ef4581)
1*f3ef73ceSBarry Smith /*$Id: mpidense.c,v 1.149 2001/01/19 23:20:28 balay 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__
12b0a32e0cSBarry Smith #define __FUNC__ "MatGetDiagonalBlock_MPIDense"
130de54da6SSatish Balay int MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
140de54da6SSatish Balay {
150de54da6SSatish Balay   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
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__
42b0a32e0cSBarry Smith #define __FUNC__ "MatSetValues_MPIDense"
438f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
448965ea79SLois Curfman McInnes {
4539b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
4639b7565bSBarry Smith   int          ierr,i,j,rstart = A->rstart,rend = A->rend,row;
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__
78b0a32e0cSBarry Smith #define __FUNC__ "MatGetValues_MPIDense"
798f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
80b49de8d1SLois Curfman McInnes {
81b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
82b49de8d1SLois Curfman McInnes   int          ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row;
83b49de8d1SLois Curfman McInnes 
843a40ed3dSBarry Smith   PetscFunctionBegin;
85b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
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__
105b0a32e0cSBarry Smith #define __FUNC__ "MatGetArray_MPIDense"
1068f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array)
107ff14e315SSatish Balay {
108ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
109ff14e315SSatish Balay   int          ierr;
110ff14e315SSatish Balay 
1113a40ed3dSBarry Smith   PetscFunctionBegin;
112ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1133a40ed3dSBarry Smith   PetscFunctionReturn(0);
114ff14e315SSatish Balay }
115ff14e315SSatish Balay 
1165615d1e5SSatish Balay #undef __FUNC__
117b0a32e0cSBarry Smith #define __FUNC__ "MatGetSubMatrix_MPIDense"
118ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
119ca3fa75bSLois Curfman McInnes {
120ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
121ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
12272d926a5SLois Curfman McInnes   int          i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols,rank;
123ca3fa75bSLois Curfman McInnes   Scalar       *av,*bv,*v = lmat->v;
124ca3fa75bSLois Curfman McInnes   Mat          newmat;
125ca3fa75bSLois Curfman McInnes 
126ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
1277eba5e9cSLois Curfman McInnes   ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr);
128ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
129ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
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__
173b0a32e0cSBarry Smith #define __FUNC__ "MatRestoreArray_MPIDense"
1748f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array)
175ff14e315SSatish Balay {
1763a40ed3dSBarry Smith   PetscFunctionBegin;
1773a40ed3dSBarry Smith   PetscFunctionReturn(0);
178ff14e315SSatish Balay }
179ff14e315SSatish Balay 
1805615d1e5SSatish Balay #undef __FUNC__
181b0a32e0cSBarry Smith #define __FUNC__ "MatAssemblyBegin_MPIDense"
1828f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
1838965ea79SLois Curfman McInnes {
18439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
1858965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
186d36fbae8SSatish Balay   int          ierr,nstash,reallocs;
1878965ea79SLois Curfman McInnes   InsertMode   addv;
1888965ea79SLois Curfman McInnes 
1893a40ed3dSBarry Smith   PetscFunctionBegin;
1908965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
191ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
1927056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
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);
199b0a32e0cSBarry Smith   PetscLogInfo(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__
204b0a32e0cSBarry Smith #define __FUNC__ "MatAssemblyEnd_MPIDense"
2058f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2068965ea79SLois Curfman McInnes {
20739ddd567SLois Curfman McInnes   Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
2087ef1d9bdSSatish Balay   int          i,n,ierr,*row,*col,flg,j,rstart,ncols;
2097ef1d9bdSSatish Balay   Scalar       *val;
210e0fa3b82SLois Curfman McInnes   InsertMode   addv=mat->insertmode;
2118965ea79SLois Curfman McInnes 
2123a40ed3dSBarry Smith   PetscFunctionBegin;
2138965ea79SLois Curfman McInnes   /*  wait on receives */
2147ef1d9bdSSatish Balay   while (1) {
2158798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2167ef1d9bdSSatish Balay     if (!flg) break;
2178965ea79SLois Curfman McInnes 
2187ef1d9bdSSatish Balay     for (i=0; i<n;) {
2197ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2207ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2217ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2227ef1d9bdSSatish Balay       else       ncols = n-i;
2237ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2247ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2257ef1d9bdSSatish Balay       i = j;
2268965ea79SLois Curfman McInnes     }
2277ef1d9bdSSatish Balay   }
2288798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2298965ea79SLois Curfman McInnes 
23039ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
23139ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2328965ea79SLois Curfman McInnes 
2336d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
23439ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2358965ea79SLois Curfman McInnes   }
2363a40ed3dSBarry Smith   PetscFunctionReturn(0);
2378965ea79SLois Curfman McInnes }
2388965ea79SLois Curfman McInnes 
2395615d1e5SSatish Balay #undef __FUNC__
240b0a32e0cSBarry Smith #define __FUNC__ "MatZeroEntries_MPIDense"
2418f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A)
2428965ea79SLois Curfman McInnes {
2433a40ed3dSBarry Smith   int          ierr;
24439ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense*)A->data;
2453a40ed3dSBarry Smith 
2463a40ed3dSBarry Smith   PetscFunctionBegin;
2473a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2483a40ed3dSBarry Smith   PetscFunctionReturn(0);
2498965ea79SLois Curfman McInnes }
2508965ea79SLois Curfman McInnes 
2515615d1e5SSatish Balay #undef __FUNC__
252b0a32e0cSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIDense"
2538f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs)
2544e220ebcSLois Curfman McInnes {
2553a40ed3dSBarry Smith   PetscFunctionBegin;
2564e220ebcSLois Curfman McInnes   *bs = 1;
2573a40ed3dSBarry Smith   PetscFunctionReturn(0);
2584e220ebcSLois Curfman McInnes }
2594e220ebcSLois Curfman McInnes 
2608965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2618965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2628965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2633501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2648965ea79SLois Curfman McInnes    routine.
2658965ea79SLois Curfman McInnes */
2665615d1e5SSatish Balay #undef __FUNC__
267b0a32e0cSBarry Smith #define __FUNC__ "MatZeroRows_MPIDense"
2688f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag)
2698965ea79SLois Curfman McInnes {
27039ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
2718965ea79SLois Curfman McInnes   int            i,ierr,N,*rows,*owners = l->rowners,size = l->size;
27235d8aa7fSBarry 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;
28035d8aa7fSBarry 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 */
287b0a32e0cSBarry Smith   ierr  = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr);
288549d3d68SSatish Balay   ierr  = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
289549d3d68SSatish Balay   procs = nprocs + size;
290b0a32e0cSBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/
2918965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
2928965ea79SLois Curfman McInnes     idx = rows[i];
29335d8aa7fSBarry Smith     found = PETSC_FALSE;
2948965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
2958965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
29635d8aa7fSBarry 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*/
304b0a32e0cSBarry Smith   ierr   = PetscMalloc(2*size*sizeof(int),&work);CHKERRQ(ierr);
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:   */
311b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr);
312b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
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   */
321b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr);
322b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
323b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr);
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 */
344b0a32e0cSBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr);
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 */
359b0a32e0cSBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr);
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);
374b0a32e0cSBarry Smith   PetscLogObjectParent(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) {
381b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
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__
392b0a32e0cSBarry Smith #define __FUNC__ "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__
406b0a32e0cSBarry Smith #define __FUNC__ "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__
420b0a32e0cSBarry Smith #define __FUNC__ "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__
436b0a32e0cSBarry Smith #define __FUNC__ "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__
451b0a32e0cSBarry Smith #define __FUNC__ "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__
474b0a32e0cSBarry Smith #define __FUNC__ "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)
483b0a32e0cSBarry Smith   PetscLogObjectState((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__
501b0a32e0cSBarry Smith #define __FUNC__ "MatView_MPIDense_Binary"
502b0a32e0cSBarry Smith static int MatView_MPIDense_Binary(Mat mat,PetscViewer 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__
516b0a32e0cSBarry Smith #define __FUNC__ "MatView_MPIDense_ASCIIorDraworSocket"
517b0a32e0cSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5188965ea79SLois Curfman McInnes {
51939ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
520fb9695e5SSatish Balay   int               ierr,size = mdn->size,rank = mdn->rank;
521b0a32e0cSBarry Smith   PetscViewerType   vtype;
522f1af5d2fSBarry Smith   PetscTruth        isascii,isdraw;
523b0a32e0cSBarry Smith   PetscViewer       sviewer;
524*f3ef73ceSBarry Smith   PetscViewerFormat format;
5258965ea79SLois Curfman McInnes 
5263a40ed3dSBarry Smith   PetscFunctionBegin;
527b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
528fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
529f1af5d2fSBarry Smith   if (isascii) {
530b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
531b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
532fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_INFO_LONG) {
5334e220ebcSLois Curfman McInnes       MatInfo info;
534888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
535b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m,
5366831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
537b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5383501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5393a40ed3dSBarry Smith       PetscFunctionReturn(0);
540fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5413a40ed3dSBarry Smith       PetscFunctionReturn(0);
5428965ea79SLois Curfman McInnes     }
543f1af5d2fSBarry Smith   } else if (isdraw) {
544b0a32e0cSBarry Smith     PetscDraw       draw;
545f1af5d2fSBarry Smith     PetscTruth isnull;
546f1af5d2fSBarry Smith 
547b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
548b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
549f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
550f1af5d2fSBarry Smith   }
55177ed5343SBarry Smith 
5528965ea79SLois Curfman McInnes   if (size == 1) {
55339ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5543a40ed3dSBarry Smith   } else {
5558965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5568965ea79SLois Curfman McInnes     Mat          A;
557273d9f13SBarry Smith     int          M = mat->M,N = mat->N,m,row,i,nz,*cols;
55839ddd567SLois Curfman McInnes     Scalar       *vals;
5598965ea79SLois Curfman McInnes 
5608965ea79SLois Curfman McInnes     if (!rank) {
561f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5623a40ed3dSBarry Smith     } else {
563f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5648965ea79SLois Curfman McInnes     }
565b0a32e0cSBarry Smith     PetscLogObjectParent(mat,A);
5668965ea79SLois Curfman McInnes 
56739ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
56839ddd567SLois Curfman McInnes        but it's quick for now */
569273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
5708965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
57139ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57239ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
57339ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57439ddd567SLois Curfman McInnes       row++;
5758965ea79SLois Curfman McInnes     }
5768965ea79SLois Curfman McInnes 
5776d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5786d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
579b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
580b9b97703SBarry Smith     if (!rank) {
5816831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5828965ea79SLois Curfman McInnes     }
583b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
584b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5858965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
5868965ea79SLois Curfman McInnes   }
5873a40ed3dSBarry Smith   PetscFunctionReturn(0);
5888965ea79SLois Curfman McInnes }
5898965ea79SLois Curfman McInnes 
5905615d1e5SSatish Balay #undef __FUNC__
591b0a32e0cSBarry Smith #define __FUNC__ "MatView_MPIDense"
592b0a32e0cSBarry Smith int MatView_MPIDense(Mat mat,PetscViewer viewer)
5938965ea79SLois Curfman McInnes {
59439ddd567SLois Curfman McInnes   int        ierr;
595f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
5968965ea79SLois Curfman McInnes 
597433994e6SBarry Smith   PetscFunctionBegin;
5980f5bd95cSBarry Smith 
599b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
600fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
601b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
602fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6030f5bd95cSBarry Smith 
604f1af5d2fSBarry Smith   if (isascii || issocket || isdraw) {
605f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6060f5bd95cSBarry Smith   } else if (isbinary) {
6073a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6085cd90555SBarry Smith   } else {
60929bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6108965ea79SLois Curfman McInnes   }
6113a40ed3dSBarry Smith   PetscFunctionReturn(0);
6128965ea79SLois Curfman McInnes }
6138965ea79SLois Curfman McInnes 
6145615d1e5SSatish Balay #undef __FUNC__
615b0a32e0cSBarry Smith #define __FUNC__ "MatGetInfo_MPIDense"
6168f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6178965ea79SLois Curfman McInnes {
6183501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6193501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6204e220ebcSLois Curfman McInnes   int          ierr;
621329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6228965ea79SLois Curfman McInnes 
6233a40ed3dSBarry Smith   PetscFunctionBegin;
624273d9f13SBarry Smith   info->rows_global    = (double)A->M;
625273d9f13SBarry Smith   info->columns_global = (double)A->N;
626273d9f13SBarry Smith   info->rows_local     = (double)A->m;
627273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6284e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6294e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6304e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6314e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6328965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6334e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6344e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6354e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6364e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6374e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6388965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
639f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
6404e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6414e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6424e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6434e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6444e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6458965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
646f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
6474e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6484e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6494e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6504e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6514e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6528965ea79SLois Curfman McInnes   }
6534e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6544e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6554e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6563a40ed3dSBarry Smith   PetscFunctionReturn(0);
6578965ea79SLois Curfman McInnes }
6588965ea79SLois Curfman McInnes 
6595615d1e5SSatish Balay #undef __FUNC__
660b0a32e0cSBarry Smith #define __FUNC__ "MatSetOption_MPIDense"
6618f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6628965ea79SLois Curfman McInnes {
66339ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
664273d9f13SBarry Smith   int          ierr;
6658965ea79SLois Curfman McInnes 
6663a40ed3dSBarry Smith   PetscFunctionBegin;
6676d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6686d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
6694787f768SSatish Balay       op == MAT_NEW_NONZERO_LOCATION_ERR ||
6704787f768SSatish Balay       op == MAT_NEW_NONZERO_ALLOCATION_ERR ||
671219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
672219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED) {
673273d9f13SBarry Smith         ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
674b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
675273d9f13SBarry Smith         a->roworiented = PETSC_TRUE;
676273d9f13SBarry Smith         ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
677b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
678219d9a1aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
6796d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
6806d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
681b51ba29fSSatish Balay              op == MAT_YES_NEW_DIAGONALS ||
682b51ba29fSSatish Balay              op == MAT_USE_HASH_TABLE) {
683b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
6843a40ed3dSBarry Smith   } else if (op == MAT_COLUMN_ORIENTED) {
685273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
686273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
6873782ba37SSatish Balay   } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) {
688273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
6893a40ed3dSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS) {
69029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
6913a40ed3dSBarry Smith   } else {
69229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
6933a40ed3dSBarry Smith   }
6943a40ed3dSBarry Smith   PetscFunctionReturn(0);
6958965ea79SLois Curfman McInnes }
6968965ea79SLois Curfman McInnes 
6975615d1e5SSatish Balay #undef __FUNC__
698b0a32e0cSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIDense"
6998f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
7008965ea79SLois Curfman McInnes {
7013501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7023a40ed3dSBarry Smith 
7033a40ed3dSBarry Smith   PetscFunctionBegin;
7044c49b128SBarry Smith   if (m) *m = mat->rstart;
7054c49b128SBarry Smith   if (n) *n = mat->rend;
7063a40ed3dSBarry Smith   PetscFunctionReturn(0);
7078965ea79SLois Curfman McInnes }
7088965ea79SLois Curfman McInnes 
7095615d1e5SSatish Balay #undef __FUNC__
710b0a32e0cSBarry Smith #define __FUNC__ "MatGetRow_MPIDense"
7118f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
7128965ea79SLois Curfman McInnes {
7133501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7143a40ed3dSBarry Smith   int          lrow,rstart = mat->rstart,rend = mat->rend,ierr;
7158965ea79SLois Curfman McInnes 
7163a40ed3dSBarry Smith   PetscFunctionBegin;
71729bbc08cSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
7188965ea79SLois Curfman McInnes   lrow = row - rstart;
7193a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7203a40ed3dSBarry Smith   PetscFunctionReturn(0);
7218965ea79SLois Curfman McInnes }
7228965ea79SLois Curfman McInnes 
7235615d1e5SSatish Balay #undef __FUNC__
724b0a32e0cSBarry Smith #define __FUNC__ "MatRestoreRow_MPIDense"
7258f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
7268965ea79SLois Curfman McInnes {
727606d414cSSatish Balay   int ierr;
728606d414cSSatish Balay 
7293a40ed3dSBarry Smith   PetscFunctionBegin;
730606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
731606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7323a40ed3dSBarry Smith   PetscFunctionReturn(0);
7338965ea79SLois Curfman McInnes }
7348965ea79SLois Curfman McInnes 
7355615d1e5SSatish Balay #undef __FUNC__
736b0a32e0cSBarry Smith #define __FUNC__ "MatDiagonalScale_MPIDense"
7375b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7385b2fa520SLois Curfman McInnes {
7395b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7405b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7415b2fa520SLois Curfman McInnes   Scalar       *l,*r,x,*v;
742273d9f13SBarry Smith   int          ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7435b2fa520SLois Curfman McInnes 
7445b2fa520SLois Curfman McInnes   PetscFunctionBegin;
74572d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7465b2fa520SLois Curfman McInnes   if (ll) {
74772d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
74829bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
7495b2fa520SLois Curfman McInnes     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7505b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7515b2fa520SLois Curfman McInnes       x = l[i];
7525b2fa520SLois Curfman McInnes       v = mat->v + i;
7535b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7545b2fa520SLois Curfman McInnes     }
7555b2fa520SLois Curfman McInnes     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
756b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7575b2fa520SLois Curfman McInnes   }
7585b2fa520SLois Curfman McInnes   if (rr) {
75972d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
76029bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7615b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7625b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7635b2fa520SLois Curfman McInnes     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
7645b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7655b2fa520SLois Curfman McInnes       x = r[i];
7665b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7675b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7685b2fa520SLois Curfman McInnes     }
76972d926a5SLois Curfman McInnes     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
770b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7715b2fa520SLois Curfman McInnes   }
7725b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7735b2fa520SLois Curfman McInnes }
7745b2fa520SLois Curfman McInnes 
7755b2fa520SLois Curfman McInnes #undef __FUNC__
776b0a32e0cSBarry Smith #define __FUNC__ "MatNorm_MPIDense"
777329f5518SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *norm)
778096963f5SLois Curfman McInnes {
7793501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7803501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7813501a2bdSLois Curfman McInnes   int          ierr,i,j;
782329f5518SBarry Smith   PetscReal    sum = 0.0;
7833501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
7843501a2bdSLois Curfman McInnes 
7853a40ed3dSBarry Smith   PetscFunctionBegin;
7863501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
7873501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm);CHKERRQ(ierr);
7883501a2bdSLois Curfman McInnes   } else {
7893501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
790273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
791aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
792329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
7933501a2bdSLois Curfman McInnes #else
7943501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7953501a2bdSLois Curfman McInnes #endif
7963501a2bdSLois Curfman McInnes       }
797ca161407SBarry Smith       ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
7983501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
799b0a32e0cSBarry Smith       PetscLogFlops(2*mdn->A->n*mdn->A->m);
8003a40ed3dSBarry Smith     } else if (type == NORM_1) {
801329f5518SBarry Smith       PetscReal *tmp,*tmp2;
802b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
803273d9f13SBarry Smith       tmp2 = tmp + A->N;
804273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
805096963f5SLois Curfman McInnes       *norm = 0.0;
8063501a2bdSLois Curfman McInnes       v = mat->v;
807273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
808273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
80967e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8103501a2bdSLois Curfman McInnes         }
8113501a2bdSLois Curfman McInnes       }
812273d9f13SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
813273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
8143501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
8153501a2bdSLois Curfman McInnes       }
816606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
817b0a32e0cSBarry Smith       PetscLogFlops(A->n*A->m);
8183a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
819329f5518SBarry Smith       PetscReal ntemp;
8203501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
821ca161407SBarry Smith       ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
8223a40ed3dSBarry Smith     } else {
82329bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8243501a2bdSLois Curfman McInnes     }
8253501a2bdSLois Curfman McInnes   }
8263a40ed3dSBarry Smith   PetscFunctionReturn(0);
8273501a2bdSLois Curfman McInnes }
8283501a2bdSLois Curfman McInnes 
8295615d1e5SSatish Balay #undef __FUNC__
830b0a32e0cSBarry Smith #define __FUNC__ "MatTranspose_MPIDense"
8318f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8323501a2bdSLois Curfman McInnes {
8333501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8343501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8353501a2bdSLois Curfman McInnes   Mat          B;
836273d9f13SBarry Smith   int          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8373501a2bdSLois Curfman McInnes   int          j,i,ierr;
8383501a2bdSLois Curfman McInnes   Scalar       *v;
8393501a2bdSLois Curfman McInnes 
8403a40ed3dSBarry Smith   PetscFunctionBegin;
8417c922b88SBarry Smith   if (!matout && M != N) {
84229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8437056b6fcSBarry Smith   }
8447056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8453501a2bdSLois Curfman McInnes 
846273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
847b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr);
8483501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8493501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8503501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8513501a2bdSLois Curfman McInnes     v   += m;
8523501a2bdSLois Curfman McInnes   }
853606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8546d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8556d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8567c922b88SBarry Smith   if (matout) {
8573501a2bdSLois Curfman McInnes     *matout = B;
8583501a2bdSLois Curfman McInnes   } else {
859273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8603501a2bdSLois Curfman McInnes   }
8613a40ed3dSBarry Smith   PetscFunctionReturn(0);
862096963f5SLois Curfman McInnes }
863096963f5SLois Curfman McInnes 
864d9eff348SSatish Balay #include "petscblaslapack.h"
8655615d1e5SSatish Balay #undef __FUNC__
866b0a32e0cSBarry Smith #define __FUNC__ "MatScale_MPIDense"
8678f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA)
86844cd7ae7SLois Curfman McInnes {
86944cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
87044cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
87144cd7ae7SLois Curfman McInnes   int          one = 1,nz;
87244cd7ae7SLois Curfman McInnes 
8733a40ed3dSBarry Smith   PetscFunctionBegin;
874273d9f13SBarry Smith   nz = inA->m*inA->N;
87544cd7ae7SLois Curfman McInnes   BLscal_(&nz,alpha,a->v,&one);
876b0a32e0cSBarry Smith   PetscLogFlops(nz);
8773a40ed3dSBarry Smith   PetscFunctionReturn(0);
87844cd7ae7SLois Curfman McInnes }
87944cd7ae7SLois Curfman McInnes 
8805609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
881ca44d042SBarry Smith EXTERN int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **);
8828965ea79SLois Curfman McInnes 
883273d9f13SBarry Smith #undef __FUNC__
884b0a32e0cSBarry Smith #define __FUNC__ "MatSetUpPreallocation_MPIDense"
885273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A)
886273d9f13SBarry Smith {
887273d9f13SBarry Smith   int        ierr;
888273d9f13SBarry Smith 
889273d9f13SBarry Smith   PetscFunctionBegin;
890273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
891273d9f13SBarry Smith   PetscFunctionReturn(0);
892273d9f13SBarry Smith }
893273d9f13SBarry Smith 
8948965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
89509dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
89609dc0095SBarry Smith        MatGetRow_MPIDense,
89709dc0095SBarry Smith        MatRestoreRow_MPIDense,
89809dc0095SBarry Smith        MatMult_MPIDense,
89909dc0095SBarry Smith        MatMultAdd_MPIDense,
9007c922b88SBarry Smith        MatMultTranspose_MPIDense,
9017c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9028965ea79SLois Curfman McInnes        0,
90309dc0095SBarry Smith        0,
90409dc0095SBarry Smith        0,
90509dc0095SBarry Smith        0,
90609dc0095SBarry Smith        0,
90709dc0095SBarry Smith        0,
90809dc0095SBarry Smith        0,
90909dc0095SBarry Smith        MatTranspose_MPIDense,
91009dc0095SBarry Smith        MatGetInfo_MPIDense,0,
91109dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9125b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
91309dc0095SBarry Smith        MatNorm_MPIDense,
91409dc0095SBarry Smith        MatAssemblyBegin_MPIDense,
91509dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
91609dc0095SBarry Smith        0,
91709dc0095SBarry Smith        MatSetOption_MPIDense,
91809dc0095SBarry Smith        MatZeroEntries_MPIDense,
91909dc0095SBarry Smith        MatZeroRows_MPIDense,
92009dc0095SBarry Smith        0,
92109dc0095SBarry Smith        0,
92209dc0095SBarry Smith        0,
92309dc0095SBarry Smith        0,
924273d9f13SBarry Smith        MatSetUpPreallocation_MPIDense,
925273d9f13SBarry Smith        0,
92639ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
92709dc0095SBarry Smith        0,
92809dc0095SBarry Smith        0,
92909dc0095SBarry Smith        MatGetArray_MPIDense,
93009dc0095SBarry Smith        MatRestoreArray_MPIDense,
9315609ef8eSBarry Smith        MatDuplicate_MPIDense,
93209dc0095SBarry Smith        0,
93309dc0095SBarry Smith        0,
93409dc0095SBarry Smith        0,
93509dc0095SBarry Smith        0,
93609dc0095SBarry Smith        0,
9372ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
93809dc0095SBarry Smith        0,
93909dc0095SBarry Smith        MatGetValues_MPIDense,
94009dc0095SBarry Smith        0,
94109dc0095SBarry Smith        0,
94209dc0095SBarry Smith        MatScale_MPIDense,
94309dc0095SBarry Smith        0,
94409dc0095SBarry Smith        0,
94509dc0095SBarry Smith        0,
94609dc0095SBarry Smith        MatGetBlockSize_MPIDense,
94709dc0095SBarry Smith        0,
94809dc0095SBarry Smith        0,
94909dc0095SBarry Smith        0,
95009dc0095SBarry Smith        0,
95109dc0095SBarry Smith        0,
95209dc0095SBarry Smith        0,
95309dc0095SBarry Smith        0,
95409dc0095SBarry Smith        0,
95509dc0095SBarry Smith        0,
956ca3fa75bSLois Curfman McInnes        MatGetSubMatrix_MPIDense,
957b9b97703SBarry Smith        MatDestroy_MPIDense,
958b9b97703SBarry Smith        MatView_MPIDense,
95909dc0095SBarry Smith        MatGetMaps_Petsc};
9608965ea79SLois Curfman McInnes 
961273d9f13SBarry Smith EXTERN_C_BEGIN
962273d9f13SBarry Smith #undef __FUNC__
963b0a32e0cSBarry Smith #define __FUNC__ "MatCreate_MPIDense"
964273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat)
965273d9f13SBarry Smith {
966273d9f13SBarry Smith   Mat_MPIDense *a;
967273d9f13SBarry Smith   int          ierr,i;
968273d9f13SBarry Smith 
969273d9f13SBarry Smith   PetscFunctionBegin;
970b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
971b0a32e0cSBarry Smith   mat->data         = (void*)a;
972273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
973273d9f13SBarry Smith   mat->factor       = 0;
974273d9f13SBarry Smith   mat->mapping      = 0;
975273d9f13SBarry Smith 
976273d9f13SBarry Smith   a->factor       = 0;
977273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
978273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
979273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
980273d9f13SBarry Smith 
981273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
982273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
983273d9f13SBarry Smith   a->nvec = mat->n;
984273d9f13SBarry Smith 
985273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
986273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
987273d9f13SBarry Smith   ierr = MapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
988273d9f13SBarry Smith   ierr = MapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
989273d9f13SBarry Smith 
990273d9f13SBarry Smith   /* build local table of row and column ownerships */
991b0a32e0cSBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr);
992273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
993b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
994273d9f13SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
995273d9f13SBarry Smith   a->rowners[0] = 0;
996273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
997273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
998273d9f13SBarry Smith   }
999273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
1000273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
1001273d9f13SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1002273d9f13SBarry Smith   a->cowners[0] = 0;
1003273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1004273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1005273d9f13SBarry Smith   }
1006273d9f13SBarry Smith 
1007273d9f13SBarry Smith   /* build cache for off array entries formed */
1008273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1009273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1010273d9f13SBarry Smith 
1011273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1012273d9f13SBarry Smith   a->lvec        = 0;
1013273d9f13SBarry Smith   a->Mvctx       = 0;
1014273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1015273d9f13SBarry Smith 
1016273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1017273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1018273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1019273d9f13SBarry Smith   PetscFunctionReturn(0);
1020273d9f13SBarry Smith }
1021273d9f13SBarry Smith EXTERN_C_END
1022273d9f13SBarry Smith 
1023273d9f13SBarry Smith #undef __FUNC__
1024b0a32e0cSBarry Smith #define __FUNC__ "MatMPIDenseSetPreallocation"
1025273d9f13SBarry Smith /*@C
1026273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1027273d9f13SBarry Smith 
1028273d9f13SBarry Smith    Not collective
1029273d9f13SBarry Smith 
1030273d9f13SBarry Smith    Input Parameters:
1031273d9f13SBarry Smith .  A - the matrix
1032273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1033273d9f13SBarry Smith    to control all matrix memory allocation.
1034273d9f13SBarry Smith 
1035273d9f13SBarry Smith    Notes:
1036273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1037273d9f13SBarry Smith    storage by columns.
1038273d9f13SBarry Smith 
1039273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1040273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1041273d9f13SBarry Smith    set data=PETSC_NULL.
1042273d9f13SBarry Smith 
1043273d9f13SBarry Smith    Level: intermediate
1044273d9f13SBarry Smith 
1045273d9f13SBarry Smith .keywords: matrix,dense, parallel
1046273d9f13SBarry Smith 
1047273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1048273d9f13SBarry Smith @*/
1049273d9f13SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,Scalar *data)
1050273d9f13SBarry Smith {
1051273d9f13SBarry Smith   Mat_MPIDense *a;
1052273d9f13SBarry Smith   int          ierr;
1053273d9f13SBarry Smith   PetscTruth   flg2;
1054273d9f13SBarry Smith 
1055273d9f13SBarry Smith   PetscFunctionBegin;
1056273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)mat,MATMPIDENSE,&flg2);CHKERRQ(ierr);
1057273d9f13SBarry Smith   if (!flg2) PetscFunctionReturn(0);
1058273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
1059273d9f13SBarry Smith   /* Note:  For now, when data is specified above, this assumes the user correctly
1060273d9f13SBarry Smith    allocates the local dense storage space.  We should add error checking. */
1061273d9f13SBarry Smith 
1062273d9f13SBarry Smith   a    = (Mat_MPIDense*)mat->data;
1063273d9f13SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mat->m,mat->N,data,&a->A);CHKERRQ(ierr);
1064b0a32e0cSBarry Smith   PetscLogObjectParent(mat,a->A);
1065273d9f13SBarry Smith   PetscFunctionReturn(0);
1066273d9f13SBarry Smith }
1067273d9f13SBarry Smith 
10685615d1e5SSatish Balay #undef __FUNC__
1069b0a32e0cSBarry Smith #define __FUNC__ "MatCreateMPIDense"
10708965ea79SLois Curfman McInnes /*@C
107139ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
10728965ea79SLois Curfman McInnes 
1073db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1074db81eaa0SLois Curfman McInnes 
10758965ea79SLois Curfman McInnes    Input Parameters:
1076db81eaa0SLois Curfman McInnes +  comm - MPI communicator
10778965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1078db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
10798965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1080db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1081db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1082dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
10838965ea79SLois Curfman McInnes 
10848965ea79SLois Curfman McInnes    Output Parameter:
1085477f1c0bSLois Curfman McInnes .  A - the matrix
10868965ea79SLois Curfman McInnes 
1087b259b22eSLois Curfman McInnes    Notes:
108839ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
108939ddd567SLois Curfman McInnes    storage by columns.
10908965ea79SLois Curfman McInnes 
109118f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
109218f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1093b4fd4287SBarry Smith    set data=PETSC_NULL.
109418f449edSLois Curfman McInnes 
10958965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
10968965ea79SLois Curfman McInnes    (possibly both).
10978965ea79SLois Curfman McInnes 
1098027ccd11SLois Curfman McInnes    Level: intermediate
1099027ccd11SLois Curfman McInnes 
110039ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
11018965ea79SLois Curfman McInnes 
110239ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
11038965ea79SLois Curfman McInnes @*/
1104477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
11058965ea79SLois Curfman McInnes {
1106273d9f13SBarry Smith   int ierr,size;
11078965ea79SLois Curfman McInnes 
11083a40ed3dSBarry Smith   PetscFunctionBegin;
1109273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1110273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1111273d9f13SBarry Smith   if (size > 1) {
1112273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1113273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1114273d9f13SBarry Smith   } else {
1115273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1116273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
11178c469469SLois Curfman McInnes   }
11183a40ed3dSBarry Smith   PetscFunctionReturn(0);
11198965ea79SLois Curfman McInnes }
11208965ea79SLois Curfman McInnes 
11215615d1e5SSatish Balay #undef __FUNC__
1122b0a32e0cSBarry Smith #define __FUNC__ "MatDuplicate_MPIDense"
11235609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11248965ea79SLois Curfman McInnes {
11258965ea79SLois Curfman McInnes   Mat          mat;
11263501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
112739ddd567SLois Curfman McInnes   int          ierr;
11288965ea79SLois Curfman McInnes 
11293a40ed3dSBarry Smith   PetscFunctionBegin;
11308965ea79SLois Curfman McInnes   *newmat       = 0;
1131273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1132273d9f13SBarry Smith   ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr);
1133b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1134b0a32e0cSBarry Smith   mat->data         = (void*)a;
1135549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
11363501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1137c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1138273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
11398965ea79SLois Curfman McInnes 
11408965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
11418965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
11428965ea79SLois Curfman McInnes   a->size         = oldmat->size;
11438965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1144e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1145b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
11463782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1147b0a32e0cSBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1148b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1149549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
11508798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
11518965ea79SLois Curfman McInnes 
1152329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
11535609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
1154b0a32e0cSBarry Smith   PetscLogObjectParent(mat,a->A);
11558965ea79SLois Curfman McInnes   *newmat = mat;
11563a40ed3dSBarry Smith   PetscFunctionReturn(0);
11578965ea79SLois Curfman McInnes }
11588965ea79SLois Curfman McInnes 
1159e090d566SSatish Balay #include "petscsys.h"
11608965ea79SLois Curfman McInnes 
11615615d1e5SSatish Balay #undef __FUNC__
1162b0a32e0cSBarry Smith #define __FUNC__ "MatLoad_MPIDense_DenseInFile"
116390ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat)
116490ace30eSBarry Smith {
116540011551SBarry Smith   int        *rowners,i,size,rank,m,ierr,nz,j;
116690ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
116790ace30eSBarry Smith   MPI_Status status;
116890ace30eSBarry Smith 
11693a40ed3dSBarry Smith   PetscFunctionBegin;
1170d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1171d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
117290ace30eSBarry Smith 
117390ace30eSBarry Smith   /* determine ownership of all rows */
117490ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
1175b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1176ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
117790ace30eSBarry Smith   rowners[0] = 0;
117890ace30eSBarry Smith   for (i=2; i<=size; i++) {
117990ace30eSBarry Smith     rowners[i] += rowners[i-1];
118090ace30eSBarry Smith   }
118190ace30eSBarry Smith 
118290ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
118390ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
118490ace30eSBarry Smith 
118590ace30eSBarry Smith   if (!rank) {
1186b0a32e0cSBarry Smith     ierr = PetscMalloc(m*N*sizeof(Scalar),&vals);CHKERRQ(ierr);
118790ace30eSBarry Smith 
118890ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
11890752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
119090ace30eSBarry Smith 
119190ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
119290ace30eSBarry Smith     vals_ptr = vals;
119390ace30eSBarry Smith     for (i=0; i<m; i++) {
119490ace30eSBarry Smith       for (j=0; j<N; j++) {
119590ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
119690ace30eSBarry Smith       }
119790ace30eSBarry Smith     }
119890ace30eSBarry Smith 
119990ace30eSBarry Smith     /* read in other processors and ship out */
120090ace30eSBarry Smith     for (i=1; i<size; i++) {
120190ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12020752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1203ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
120490ace30eSBarry Smith     }
12053a40ed3dSBarry Smith   } else {
120690ace30eSBarry Smith     /* receive numeric values */
1207b0a32e0cSBarry Smith     ierr = PetscMalloc(m*N*sizeof(Scalar),&vals);CHKERRQ(ierr);
120890ace30eSBarry Smith 
120990ace30eSBarry Smith     /* receive message of values*/
1210ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
121190ace30eSBarry Smith 
121290ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
121390ace30eSBarry Smith     vals_ptr = vals;
121490ace30eSBarry Smith     for (i=0; i<m; i++) {
121590ace30eSBarry Smith       for (j=0; j<N; j++) {
121690ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
121790ace30eSBarry Smith       }
121890ace30eSBarry Smith     }
121990ace30eSBarry Smith   }
1220606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1221606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
12226d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12236d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12243a40ed3dSBarry Smith   PetscFunctionReturn(0);
122590ace30eSBarry Smith }
122690ace30eSBarry Smith 
1227273d9f13SBarry Smith EXTERN_C_BEGIN
12285615d1e5SSatish Balay #undef __FUNC__
1229b0a32e0cSBarry Smith #define __FUNC__ "MatLoad_MPIDense"
1230b0a32e0cSBarry Smith int MatLoad_MPIDense(PetscViewer viewer,MatType type,Mat *newmat)
12318965ea79SLois Curfman McInnes {
12328965ea79SLois Curfman McInnes   Mat          A;
12338965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
123419bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
12358965ea79SLois Curfman McInnes   MPI_Status   status;
12368965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
12378965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols;
123819bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
12393a40ed3dSBarry Smith   int          i,nz,ierr,j,rstart,rend,fd;
12408965ea79SLois Curfman McInnes 
12413a40ed3dSBarry Smith   PetscFunctionBegin;
1242d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1243d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12448965ea79SLois Curfman McInnes   if (!rank) {
1245b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
12460752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
124729bbc08cSBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
12488965ea79SLois Curfman McInnes   }
12498965ea79SLois Curfman McInnes 
1250ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
125190ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
125290ace30eSBarry Smith 
125390ace30eSBarry Smith   /*
125490ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
125590ace30eSBarry Smith   */
125690ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
12573a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
12583a40ed3dSBarry Smith     PetscFunctionReturn(0);
125990ace30eSBarry Smith   }
126090ace30eSBarry Smith 
12618965ea79SLois Curfman McInnes   /* determine ownership of all rows */
12628965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
1263b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1264ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
12658965ea79SLois Curfman McInnes   rowners[0] = 0;
12668965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
12678965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
12688965ea79SLois Curfman McInnes   }
12698965ea79SLois Curfman McInnes   rstart = rowners[rank];
12708965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
12718965ea79SLois Curfman McInnes 
12728965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
1273b0a32e0cSBarry Smith   ierr    = PetscMalloc(2*(rend-rstart)*sizeof(int),&ourlens);CHKERRQ(ierr);
12748965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
12758965ea79SLois Curfman McInnes   if (!rank) {
1276b0a32e0cSBarry Smith     ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr);
12770752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1278b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr);
12798965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1280ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1281606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1282ca161407SBarry Smith   } else {
1283ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
12848965ea79SLois Curfman McInnes   }
12858965ea79SLois Curfman McInnes 
12868965ea79SLois Curfman McInnes   if (!rank) {
12878965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
1288b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr);
1289549d3d68SSatish Balay     ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
12908965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
12918965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
12928965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
12938965ea79SLois Curfman McInnes       }
12948965ea79SLois Curfman McInnes     }
1295606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
12968965ea79SLois Curfman McInnes 
12978965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
12988965ea79SLois Curfman McInnes     maxnz = 0;
12998965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13000452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
13018965ea79SLois Curfman McInnes     }
1302b0a32e0cSBarry Smith     ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr);
13038965ea79SLois Curfman McInnes 
13048965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13058965ea79SLois Curfman McInnes     nz = procsnz[0];
1306b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13070752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13088965ea79SLois Curfman McInnes 
13098965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13108965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13118965ea79SLois Curfman McInnes       nz   = procsnz[i];
13120752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1313ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13148965ea79SLois Curfman McInnes     }
1315606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13163a40ed3dSBarry Smith   } else {
13178965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13188965ea79SLois Curfman McInnes     nz = 0;
13198965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13208965ea79SLois Curfman McInnes       nz += ourlens[i];
13218965ea79SLois Curfman McInnes     }
1322b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13238965ea79SLois Curfman McInnes 
13248965ea79SLois Curfman McInnes     /* receive message of column indices*/
1325ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1326ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
132729bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13288965ea79SLois Curfman McInnes   }
13298965ea79SLois Curfman McInnes 
13308965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1331549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
13328965ea79SLois Curfman McInnes   jj = 0;
13338965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13348965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
13358965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
13368965ea79SLois Curfman McInnes       jj++;
13378965ea79SLois Curfman McInnes     }
13388965ea79SLois Curfman McInnes   }
13398965ea79SLois Curfman McInnes 
13408965ea79SLois Curfman McInnes   /* create our matrix */
13418965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13428965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
13438965ea79SLois Curfman McInnes   }
1344b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
13458965ea79SLois Curfman McInnes   A = *newmat;
13468965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
13478965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
13488965ea79SLois Curfman McInnes   }
13498965ea79SLois Curfman McInnes 
13508965ea79SLois Curfman McInnes   if (!rank) {
1351b0a32e0cSBarry Smith     ierr = PetscMalloc(maxnz*sizeof(Scalar),&vals);CHKERRQ(ierr);
13528965ea79SLois Curfman McInnes 
13538965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
13548965ea79SLois Curfman McInnes     nz = procsnz[0];
13550752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
13568965ea79SLois Curfman McInnes 
13578965ea79SLois Curfman McInnes     /* insert into matrix */
13588965ea79SLois Curfman McInnes     jj      = rstart;
13598965ea79SLois Curfman McInnes     smycols = mycols;
13608965ea79SLois Curfman McInnes     svals   = vals;
13618965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13628965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13638965ea79SLois Curfman McInnes       smycols += ourlens[i];
13648965ea79SLois Curfman McInnes       svals   += ourlens[i];
13658965ea79SLois Curfman McInnes       jj++;
13668965ea79SLois Curfman McInnes     }
13678965ea79SLois Curfman McInnes 
13688965ea79SLois Curfman McInnes     /* read in other processors and ship out */
13698965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13708965ea79SLois Curfman McInnes       nz   = procsnz[i];
13710752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1372ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
13738965ea79SLois Curfman McInnes     }
1374606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
13753a40ed3dSBarry Smith   } else {
13768965ea79SLois Curfman McInnes     /* receive numeric values */
1377b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(Scalar),&vals);CHKERRQ(ierr);
13788965ea79SLois Curfman McInnes 
13798965ea79SLois Curfman McInnes     /* receive message of values*/
1380ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1381ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
138229bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
13838965ea79SLois Curfman McInnes 
13848965ea79SLois Curfman McInnes     /* insert into matrix */
13858965ea79SLois Curfman McInnes     jj      = rstart;
13868965ea79SLois Curfman McInnes     smycols = mycols;
13878965ea79SLois Curfman McInnes     svals   = vals;
13888965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13898965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13908965ea79SLois Curfman McInnes       smycols += ourlens[i];
13918965ea79SLois Curfman McInnes       svals   += ourlens[i];
13928965ea79SLois Curfman McInnes       jj++;
13938965ea79SLois Curfman McInnes     }
13948965ea79SLois Curfman McInnes   }
1395606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1396606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1397606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1398606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
13998965ea79SLois Curfman McInnes 
14006d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14016d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14023a40ed3dSBarry Smith   PetscFunctionReturn(0);
14038965ea79SLois Curfman McInnes }
1404273d9f13SBarry Smith EXTERN_C_END
140590ace30eSBarry Smith 
140690ace30eSBarry Smith 
140790ace30eSBarry Smith 
140890ace30eSBarry Smith 
140990ace30eSBarry Smith 
1410