xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 878740d94d958d3c7a5a5af3e703f834afdabafb)
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);
275c5985e7SKris Buschelman   ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr);
285c5985e7SKris Buschelman   ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr);
295c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);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 
394a2ae208SSatish Balay #undef __FUNCT__
404a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
41f15d580aSBarry Smith int MatSetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],const PetscScalar v[],InsertMode addv)
428965ea79SLois Curfman McInnes {
4339b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
4439b7565bSBarry Smith   int          ierr,i,j,rstart = A->rstart,rend = A->rend,row;
45273d9f13SBarry Smith   PetscTruth   roworiented = A->roworiented;
468965ea79SLois Curfman McInnes 
473a40ed3dSBarry Smith   PetscFunctionBegin;
488965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
495ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
50273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
518965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
528965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
5339b7565bSBarry Smith       if (roworiented) {
5439b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
553a40ed3dSBarry Smith       } else {
568965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
575ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
58273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
5939b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
6039b7565bSBarry Smith         }
618965ea79SLois Curfman McInnes       }
623a40ed3dSBarry Smith     } else {
633782ba37SSatish Balay       if (!A->donotstash) {
6439b7565bSBarry Smith         if (roworiented) {
658798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
66d36fbae8SSatish Balay         } else {
678798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
6839b7565bSBarry Smith         }
69b49de8d1SLois Curfman McInnes       }
70b49de8d1SLois Curfman McInnes     }
713782ba37SSatish Balay   }
723a40ed3dSBarry Smith   PetscFunctionReturn(0);
73b49de8d1SLois Curfman McInnes }
74b49de8d1SLois Curfman McInnes 
754a2ae208SSatish Balay #undef __FUNCT__
764a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
77f15d580aSBarry Smith int MatGetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],PetscScalar v[])
78b49de8d1SLois Curfman McInnes {
79b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
80b49de8d1SLois Curfman McInnes   int          ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row;
81b49de8d1SLois Curfman McInnes 
823a40ed3dSBarry Smith   PetscFunctionBegin;
83b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
8429bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
85273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
86b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
87b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
88b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
8929bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
90273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
9129bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
92a8c6a408SBarry Smith         }
93b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
94b49de8d1SLois Curfman McInnes       }
95a8c6a408SBarry Smith     } else {
9629bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
978965ea79SLois Curfman McInnes     }
988965ea79SLois Curfman McInnes   }
993a40ed3dSBarry Smith   PetscFunctionReturn(0);
1008965ea79SLois Curfman McInnes }
1018965ea79SLois Curfman McInnes 
1024a2ae208SSatish Balay #undef __FUNCT__
1034a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
1044e7234bfSBarry Smith int MatGetArray_MPIDense(Mat A,PetscScalar *array[])
105ff14e315SSatish Balay {
106ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
107ff14e315SSatish Balay   int          ierr;
108ff14e315SSatish Balay 
1093a40ed3dSBarry Smith   PetscFunctionBegin;
110ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1113a40ed3dSBarry Smith   PetscFunctionReturn(0);
112ff14e315SSatish Balay }
113ff14e315SSatish Balay 
1144a2ae208SSatish Balay #undef __FUNCT__
1154a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
116ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
117ca3fa75bSLois Curfman McInnes {
118ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
119ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
120cfce73b9SSatish Balay   int          i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
12187828ca2SBarry Smith   PetscScalar  *av,*bv,*v = lmat->v;
122ca3fa75bSLois Curfman McInnes   Mat          newmat;
123ca3fa75bSLois Curfman McInnes 
124ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
125ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
126ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
127b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
128b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
129ca3fa75bSLois Curfman McInnes 
130ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1317eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1327eba5e9cSLois Curfman McInnes      original matrix! */
133ca3fa75bSLois Curfman McInnes 
134ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1357eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
136ca3fa75bSLois Curfman McInnes 
137ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
138ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
13929bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
1407eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
141ca3fa75bSLois Curfman McInnes     newmat = *B;
142ca3fa75bSLois Curfman McInnes   } else {
143ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
144*878740d9SKris Buschelman     ierr = MatCreate(A->comm,nrows,cs,PETSC_DECIDE,ncols,&newmat);CHKERRQ(ierr);
145*878740d9SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
146*878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
147ca3fa75bSLois Curfman McInnes   }
148ca3fa75bSLois Curfman McInnes 
149ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
150ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
151ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
152ca3fa75bSLois Curfman McInnes 
153ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
154ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
155ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
1567eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
157ca3fa75bSLois Curfman McInnes     }
158ca3fa75bSLois Curfman McInnes   }
159ca3fa75bSLois Curfman McInnes 
160ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
161ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
163ca3fa75bSLois Curfman McInnes 
164ca3fa75bSLois Curfman McInnes   /* Free work space */
165ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
166ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
167ca3fa75bSLois Curfman McInnes   *B = newmat;
168ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
169ca3fa75bSLois Curfman McInnes }
170ca3fa75bSLois Curfman McInnes 
1714a2ae208SSatish Balay #undef __FUNCT__
1724a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
1734e7234bfSBarry Smith int MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
174ff14e315SSatish Balay {
1753a40ed3dSBarry Smith   PetscFunctionBegin;
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177ff14e315SSatish Balay }
178ff14e315SSatish Balay 
1794a2ae208SSatish Balay #undef __FUNCT__
1804a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
1818f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
1828965ea79SLois Curfman McInnes {
18339ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
1848965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
185d36fbae8SSatish Balay   int          ierr,nstash,reallocs;
1868965ea79SLois Curfman McInnes   InsertMode   addv;
1878965ea79SLois Curfman McInnes 
1883a40ed3dSBarry Smith   PetscFunctionBegin;
1898965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
190ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
1917056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
19229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
1938965ea79SLois Curfman McInnes   }
194e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
1958965ea79SLois Curfman McInnes 
1968798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
1978798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
198b0a32e0cSBarry Smith   PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs);
1993a40ed3dSBarry Smith   PetscFunctionReturn(0);
2008965ea79SLois Curfman McInnes }
2018965ea79SLois Curfman McInnes 
2024a2ae208SSatish Balay #undef __FUNCT__
2034a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
2048f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2058965ea79SLois Curfman McInnes {
20639ddd567SLois Curfman McInnes   Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
2077ef1d9bdSSatish Balay   int          i,n,ierr,*row,*col,flg,j,rstart,ncols;
20887828ca2SBarry Smith   PetscScalar  *val;
209e0fa3b82SLois Curfman McInnes   InsertMode   addv=mat->insertmode;
2108965ea79SLois Curfman McInnes 
2113a40ed3dSBarry Smith   PetscFunctionBegin;
2128965ea79SLois Curfman McInnes   /*  wait on receives */
2137ef1d9bdSSatish Balay   while (1) {
2148798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2157ef1d9bdSSatish Balay     if (!flg) break;
2168965ea79SLois Curfman McInnes 
2177ef1d9bdSSatish Balay     for (i=0; i<n;) {
2187ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2197ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2207ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2217ef1d9bdSSatish Balay       else       ncols = n-i;
2227ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2237ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2247ef1d9bdSSatish Balay       i = j;
2258965ea79SLois Curfman McInnes     }
2267ef1d9bdSSatish Balay   }
2278798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2288965ea79SLois Curfman McInnes 
22939ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
23039ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2318965ea79SLois Curfman McInnes 
2326d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
23339ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2348965ea79SLois Curfman McInnes   }
2353a40ed3dSBarry Smith   PetscFunctionReturn(0);
2368965ea79SLois Curfman McInnes }
2378965ea79SLois Curfman McInnes 
2384a2ae208SSatish Balay #undef __FUNCT__
2394a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
2408f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A)
2418965ea79SLois Curfman McInnes {
2423a40ed3dSBarry Smith   int          ierr;
24339ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense*)A->data;
2443a40ed3dSBarry Smith 
2453a40ed3dSBarry Smith   PetscFunctionBegin;
2463a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2473a40ed3dSBarry Smith   PetscFunctionReturn(0);
2488965ea79SLois Curfman McInnes }
2498965ea79SLois Curfman McInnes 
2504a2ae208SSatish Balay #undef __FUNCT__
2514a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIDense"
2528f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs)
2534e220ebcSLois Curfman McInnes {
2543a40ed3dSBarry Smith   PetscFunctionBegin;
2554e220ebcSLois Curfman McInnes   *bs = 1;
2563a40ed3dSBarry Smith   PetscFunctionReturn(0);
2574e220ebcSLois Curfman McInnes }
2584e220ebcSLois Curfman McInnes 
2598965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2608965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2618965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2623501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2638965ea79SLois Curfman McInnes    routine.
2648965ea79SLois Curfman McInnes */
2654a2ae208SSatish Balay #undef __FUNCT__
2664a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
267268466fbSBarry Smith int MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag)
2688965ea79SLois Curfman McInnes {
26939ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
2708965ea79SLois Curfman McInnes   int            i,ierr,N,*rows,*owners = l->rowners,size = l->size;
271c1dc657dSBarry Smith   int            *nprocs,j,idx,nsends;
2728965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2738965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2748965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2758965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2768965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2778965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2788965ea79SLois Curfman McInnes   IS             istmp;
27935d8aa7fSBarry Smith   PetscTruth     found;
2808965ea79SLois Curfman McInnes 
2813a40ed3dSBarry Smith   PetscFunctionBegin;
282b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
2838965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
2848965ea79SLois Curfman McInnes 
2858965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
286b0a32e0cSBarry Smith   ierr  = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr);
287549d3d68SSatish Balay   ierr  = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
288b0a32e0cSBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/
2898965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
2908965ea79SLois Curfman McInnes     idx = rows[i];
29135d8aa7fSBarry Smith     found = PETSC_FALSE;
2928965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
2938965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
294c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
2958965ea79SLois Curfman McInnes       }
2968965ea79SLois Curfman McInnes     }
29729bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
2988965ea79SLois Curfman McInnes   }
299c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3008965ea79SLois Curfman McInnes 
3018965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
302c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3038965ea79SLois Curfman McInnes 
3048965ea79SLois Curfman McInnes   /* post receives:   */
305b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr);
306b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3078965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
308ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3098965ea79SLois Curfman McInnes   }
3108965ea79SLois Curfman McInnes 
3118965ea79SLois Curfman McInnes   /* do sends:
3128965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3138965ea79SLois Curfman McInnes          the ith processor
3148965ea79SLois Curfman McInnes   */
315b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr);
316b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
317b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr);
3188965ea79SLois Curfman McInnes   starts[0]  = 0;
319c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3208965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3218965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3228965ea79SLois Curfman McInnes   }
3238965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3248965ea79SLois Curfman McInnes 
3258965ea79SLois Curfman McInnes   starts[0] = 0;
326c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3278965ea79SLois Curfman McInnes   count = 0;
3288965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
329c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
330c1dc657dSBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3318965ea79SLois Curfman McInnes     }
3328965ea79SLois Curfman McInnes   }
333606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3348965ea79SLois Curfman McInnes 
3358965ea79SLois Curfman McInnes   base = owners[rank];
3368965ea79SLois Curfman McInnes 
3378965ea79SLois Curfman McInnes   /*  wait on receives */
338b0a32e0cSBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr);
3398965ea79SLois Curfman McInnes   source = lens + nrecvs;
3408965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3418965ea79SLois Curfman McInnes   while (count) {
342ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3438965ea79SLois Curfman McInnes     /* unpack receives into our local space */
344ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3458965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3468965ea79SLois Curfman McInnes     lens[imdex]  = n;
3478965ea79SLois Curfman McInnes     slen += n;
3488965ea79SLois Curfman McInnes     count--;
3498965ea79SLois Curfman McInnes   }
350606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3518965ea79SLois Curfman McInnes 
3528965ea79SLois Curfman McInnes   /* move the data into the send scatter */
353b0a32e0cSBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr);
3548965ea79SLois Curfman McInnes   count = 0;
3558965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
3568965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3578965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
3588965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3598965ea79SLois Curfman McInnes     }
3608965ea79SLois Curfman McInnes   }
361606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
362606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
363606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
364606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3658965ea79SLois Curfman McInnes 
3668965ea79SLois Curfman McInnes   /* actually zap the local rows */
367029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
368b0a32e0cSBarry Smith   PetscLogObjectParent(A,istmp);
369606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3708965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3718965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3728965ea79SLois Curfman McInnes 
3738965ea79SLois Curfman McInnes   /* wait on sends */
3748965ea79SLois Curfman McInnes   if (nsends) {
375b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
376ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
377606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
3788965ea79SLois Curfman McInnes   }
379606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
380606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
3818965ea79SLois Curfman McInnes 
3823a40ed3dSBarry Smith   PetscFunctionReturn(0);
3838965ea79SLois Curfman McInnes }
3848965ea79SLois Curfman McInnes 
3854a2ae208SSatish Balay #undef __FUNCT__
3864a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
3878f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3888965ea79SLois Curfman McInnes {
38939ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
3908965ea79SLois Curfman McInnes   int          ierr;
391c456f294SBarry Smith 
3923a40ed3dSBarry Smith   PetscFunctionBegin;
39343a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39443a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39544cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
3963a40ed3dSBarry Smith   PetscFunctionReturn(0);
3978965ea79SLois Curfman McInnes }
3988965ea79SLois Curfman McInnes 
3994a2ae208SSatish Balay #undef __FUNCT__
4004a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
4018f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4028965ea79SLois Curfman McInnes {
40339ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4048965ea79SLois Curfman McInnes   int          ierr;
405c456f294SBarry Smith 
4063a40ed3dSBarry Smith   PetscFunctionBegin;
40743a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40843a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40944cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4103a40ed3dSBarry Smith   PetscFunctionReturn(0);
4118965ea79SLois Curfman McInnes }
4128965ea79SLois Curfman McInnes 
4134a2ae208SSatish Balay #undef __FUNCT__
4144a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
4157c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
416096963f5SLois Curfman McInnes {
417096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
418096963f5SLois Curfman McInnes   int          ierr;
41987828ca2SBarry Smith   PetscScalar  zero = 0.0;
420096963f5SLois Curfman McInnes 
4213a40ed3dSBarry Smith   PetscFunctionBegin;
4223501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
4237c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
424537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
425537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4263a40ed3dSBarry Smith   PetscFunctionReturn(0);
427096963f5SLois Curfman McInnes }
428096963f5SLois Curfman McInnes 
4294a2ae208SSatish Balay #undef __FUNCT__
4304a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
4317c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
432096963f5SLois Curfman McInnes {
433096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
434096963f5SLois Curfman McInnes   int          ierr;
435096963f5SLois Curfman McInnes 
4363a40ed3dSBarry Smith   PetscFunctionBegin;
4373501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4387c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
439537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
440537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4413a40ed3dSBarry Smith   PetscFunctionReturn(0);
442096963f5SLois Curfman McInnes }
443096963f5SLois Curfman McInnes 
4444a2ae208SSatish Balay #undef __FUNCT__
4454a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
4468f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v)
4478965ea79SLois Curfman McInnes {
44839ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
449096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
450273d9f13SBarry Smith   int          ierr,len,i,n,m = A->m,radd;
45187828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
452ed3cc1f0SBarry Smith 
4533a40ed3dSBarry Smith   PetscFunctionBegin;
454273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
455b1d4fb26SBarry Smith   ierr = VecGetArrayFast(v,&x);CHKERRQ(ierr);
456096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
457273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
458273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
4597ddc982cSLois Curfman McInnes   radd = a->rstart*m;
46044cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
461096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
462096963f5SLois Curfman McInnes   }
463b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(v,&x);CHKERRQ(ierr);
4643a40ed3dSBarry Smith   PetscFunctionReturn(0);
4658965ea79SLois Curfman McInnes }
4668965ea79SLois Curfman McInnes 
4674a2ae208SSatish Balay #undef __FUNCT__
4684a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
469e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat)
4708965ea79SLois Curfman McInnes {
4713501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4728965ea79SLois Curfman McInnes   int          ierr;
473ed3cc1f0SBarry Smith 
4743a40ed3dSBarry Smith   PetscFunctionBegin;
47594d884c6SBarry Smith 
476aa482453SBarry Smith #if defined(PETSC_USE_LOG)
477b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N);
4788965ea79SLois Curfman McInnes #endif
4798798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
480606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
4813501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
4823501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4833501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
484622d7880SLois Curfman McInnes   if (mdn->factor) {
485606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
486606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
487606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
488606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
489622d7880SLois Curfman McInnes   }
490606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
4913a40ed3dSBarry Smith   PetscFunctionReturn(0);
4928965ea79SLois Curfman McInnes }
49339ddd567SLois Curfman McInnes 
4944a2ae208SSatish Balay #undef __FUNCT__
4954a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
496b0a32e0cSBarry Smith static int MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
4978965ea79SLois Curfman McInnes {
49839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4998965ea79SLois Curfman McInnes   int          ierr;
5007056b6fcSBarry Smith 
5013a40ed3dSBarry Smith   PetscFunctionBegin;
50239ddd567SLois Curfman McInnes   if (mdn->size == 1) {
50339ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5048965ea79SLois Curfman McInnes   }
50529bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5063a40ed3dSBarry Smith   PetscFunctionReturn(0);
5078965ea79SLois Curfman McInnes }
5088965ea79SLois Curfman McInnes 
5094a2ae208SSatish Balay #undef __FUNCT__
5104a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
511b0a32e0cSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5128965ea79SLois Curfman McInnes {
51339ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
514fb9695e5SSatish Balay   int               ierr,size = mdn->size,rank = mdn->rank;
515b0a32e0cSBarry Smith   PetscViewerType   vtype;
516f1af5d2fSBarry Smith   PetscTruth        isascii,isdraw;
517b0a32e0cSBarry Smith   PetscViewer       sviewer;
518f3ef73ceSBarry Smith   PetscViewerFormat format;
5198965ea79SLois Curfman McInnes 
5203a40ed3dSBarry Smith   PetscFunctionBegin;
521b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
522fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
523f1af5d2fSBarry Smith   if (isascii) {
524b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
525b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
526456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
5274e220ebcSLois Curfman McInnes       MatInfo info;
528888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
529b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m,
5306831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
531b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5323501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5333a40ed3dSBarry Smith       PetscFunctionReturn(0);
534fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5353a40ed3dSBarry Smith       PetscFunctionReturn(0);
5368965ea79SLois Curfman McInnes     }
537f1af5d2fSBarry Smith   } else if (isdraw) {
538b0a32e0cSBarry Smith     PetscDraw       draw;
539f1af5d2fSBarry Smith     PetscTruth isnull;
540f1af5d2fSBarry Smith 
541b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
542b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
543f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
544f1af5d2fSBarry Smith   }
54577ed5343SBarry Smith 
5468965ea79SLois Curfman McInnes   if (size == 1) {
54739ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5483a40ed3dSBarry Smith   } else {
5498965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5508965ea79SLois Curfman McInnes     Mat          A;
551273d9f13SBarry Smith     int          M = mat->M,N = mat->N,m,row,i,nz,*cols;
55287828ca2SBarry Smith     PetscScalar  *vals;
5538965ea79SLois Curfman McInnes 
5548965ea79SLois Curfman McInnes     if (!rank) {
555*878740d9SKris Buschelman       ierr = MatCreate(mat->comm,M,N,M,N,&A);CHKERRQ(ierr);
5563a40ed3dSBarry Smith     } else {
557*878740d9SKris Buschelman       ierr = MatCreate(mat->comm,0,0,M,N,&A);CHKERRQ(ierr);
5588965ea79SLois Curfman McInnes     }
559*878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
560*878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
561b0a32e0cSBarry Smith     PetscLogObjectParent(mat,A);
5628965ea79SLois Curfman McInnes 
56339ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
56439ddd567SLois Curfman McInnes        but it's quick for now */
565273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
5668965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
56739ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
56839ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
56939ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57039ddd567SLois Curfman McInnes       row++;
5718965ea79SLois Curfman McInnes     }
5728965ea79SLois Curfman McInnes 
5736d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5746d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
575b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
576b9b97703SBarry Smith     if (!rank) {
5776831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5788965ea79SLois Curfman McInnes     }
579b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
580b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5818965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
5828965ea79SLois Curfman McInnes   }
5833a40ed3dSBarry Smith   PetscFunctionReturn(0);
5848965ea79SLois Curfman McInnes }
5858965ea79SLois Curfman McInnes 
5864a2ae208SSatish Balay #undef __FUNCT__
5874a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
588b0a32e0cSBarry Smith int MatView_MPIDense(Mat mat,PetscViewer viewer)
5898965ea79SLois Curfman McInnes {
59039ddd567SLois Curfman McInnes   int        ierr;
591f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
5928965ea79SLois Curfman McInnes 
593433994e6SBarry Smith   PetscFunctionBegin;
5940f5bd95cSBarry Smith 
595b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
596fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
597b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
598fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
5990f5bd95cSBarry Smith 
600f1af5d2fSBarry Smith   if (isascii || issocket || isdraw) {
601f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6020f5bd95cSBarry Smith   } else if (isbinary) {
6033a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6045cd90555SBarry Smith   } else {
60529bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6068965ea79SLois Curfman McInnes   }
6073a40ed3dSBarry Smith   PetscFunctionReturn(0);
6088965ea79SLois Curfman McInnes }
6098965ea79SLois Curfman McInnes 
6104a2ae208SSatish Balay #undef __FUNCT__
6114a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
6128f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6138965ea79SLois Curfman McInnes {
6143501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6153501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6164e220ebcSLois Curfman McInnes   int          ierr;
617329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6188965ea79SLois Curfman McInnes 
6193a40ed3dSBarry Smith   PetscFunctionBegin;
620273d9f13SBarry Smith   info->rows_global    = (double)A->M;
621273d9f13SBarry Smith   info->columns_global = (double)A->N;
622273d9f13SBarry Smith   info->rows_local     = (double)A->m;
623273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6244e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6254e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6264e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6274e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6288965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6294e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6304e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6314e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6324e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6334e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6348965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
635d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
6364e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6374e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6384e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6394e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6404e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6418965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
642d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
6434e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6444e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6454e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6464e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6474e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6488965ea79SLois Curfman McInnes   }
6494e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6504e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6514e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6523a40ed3dSBarry Smith   PetscFunctionReturn(0);
6538965ea79SLois Curfman McInnes }
6548965ea79SLois Curfman McInnes 
6554a2ae208SSatish Balay #undef __FUNCT__
6564a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
6578f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6588965ea79SLois Curfman McInnes {
65939ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
660273d9f13SBarry Smith   int          ierr;
6618965ea79SLois Curfman McInnes 
6623a40ed3dSBarry Smith   PetscFunctionBegin;
66312c028f9SKris Buschelman   switch (op) {
66412c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
66512c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
66612c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
66712c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
66812c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
66912c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
670273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
67112c028f9SKris Buschelman     break;
67212c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
673273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
674273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
67512c028f9SKris Buschelman     break;
67612c028f9SKris Buschelman   case MAT_ROWS_SORTED:
67712c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
67812c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
67912c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
680b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
68112c028f9SKris Buschelman     break;
68212c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
683273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
684273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
68512c028f9SKris Buschelman     break;
68612c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
687273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
68812c028f9SKris Buschelman     break;
68912c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
69029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
69177e54ba9SKris Buschelman   case MAT_SYMMETRIC:
69277e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
6939a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
6949a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
6959a4540c5SBarry Smith   case MAT_HERMITIAN:
6969a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
6979a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
6989a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
69977e54ba9SKris Buschelman     break;
70012c028f9SKris Buschelman   default:
70129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
7023a40ed3dSBarry Smith   }
7033a40ed3dSBarry Smith   PetscFunctionReturn(0);
7048965ea79SLois Curfman McInnes }
7058965ea79SLois Curfman McInnes 
7064a2ae208SSatish Balay #undef __FUNCT__
7074a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIDense"
70887828ca2SBarry Smith int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v)
7098965ea79SLois Curfman McInnes {
7103501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7113a40ed3dSBarry Smith   int          lrow,rstart = mat->rstart,rend = mat->rend,ierr;
7128965ea79SLois Curfman McInnes 
7133a40ed3dSBarry Smith   PetscFunctionBegin;
71429bbc08cSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
7158965ea79SLois Curfman McInnes   lrow = row - rstart;
7163a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7173a40ed3dSBarry Smith   PetscFunctionReturn(0);
7188965ea79SLois Curfman McInnes }
7198965ea79SLois Curfman McInnes 
7204a2ae208SSatish Balay #undef __FUNCT__
7214a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIDense"
72287828ca2SBarry Smith int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v)
7238965ea79SLois Curfman McInnes {
724606d414cSSatish Balay   int ierr;
725606d414cSSatish Balay 
7263a40ed3dSBarry Smith   PetscFunctionBegin;
727606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
728606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7293a40ed3dSBarry Smith   PetscFunctionReturn(0);
7308965ea79SLois Curfman McInnes }
7318965ea79SLois Curfman McInnes 
7324a2ae208SSatish Balay #undef __FUNCT__
7334a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
7345b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7355b2fa520SLois Curfman McInnes {
7365b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7375b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
73887828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
739273d9f13SBarry Smith   int          ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7405b2fa520SLois Curfman McInnes 
7415b2fa520SLois Curfman McInnes   PetscFunctionBegin;
74272d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7435b2fa520SLois Curfman McInnes   if (ll) {
74472d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
74529bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
746b1d4fb26SBarry Smith     ierr = VecGetArrayFast(ll,&l);CHKERRQ(ierr);
7475b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7485b2fa520SLois Curfman McInnes       x = l[i];
7495b2fa520SLois Curfman McInnes       v = mat->v + i;
7505b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7515b2fa520SLois Curfman McInnes     }
752b1d4fb26SBarry Smith     ierr = VecRestoreArrayFast(ll,&l);CHKERRQ(ierr);
753b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7545b2fa520SLois Curfman McInnes   }
7555b2fa520SLois Curfman McInnes   if (rr) {
75672d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
75729bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7585b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7595b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
760b1d4fb26SBarry Smith     ierr = VecGetArrayFast(mdn->lvec,&r);CHKERRQ(ierr);
7615b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7625b2fa520SLois Curfman McInnes       x = r[i];
7635b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7645b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7655b2fa520SLois Curfman McInnes     }
766b1d4fb26SBarry Smith     ierr = VecRestoreArrayFast(mdn->lvec,&r);CHKERRQ(ierr);
767b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7685b2fa520SLois Curfman McInnes   }
7695b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7705b2fa520SLois Curfman McInnes }
7715b2fa520SLois Curfman McInnes 
7724a2ae208SSatish Balay #undef __FUNCT__
7734a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
774064f8208SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
775096963f5SLois Curfman McInnes {
7763501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7773501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7783501a2bdSLois Curfman McInnes   int          ierr,i,j;
779329f5518SBarry Smith   PetscReal    sum = 0.0;
78087828ca2SBarry Smith   PetscScalar  *v = mat->v;
7813501a2bdSLois Curfman McInnes 
7823a40ed3dSBarry Smith   PetscFunctionBegin;
7833501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
784064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
7853501a2bdSLois Curfman McInnes   } else {
7863501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
787273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
788aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
789329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
7903501a2bdSLois Curfman McInnes #else
7913501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7923501a2bdSLois Curfman McInnes #endif
7933501a2bdSLois Curfman McInnes       }
794064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
795064f8208SBarry Smith       *nrm = sqrt(*nrm);
796b0a32e0cSBarry Smith       PetscLogFlops(2*mdn->A->n*mdn->A->m);
7973a40ed3dSBarry Smith     } else if (type == NORM_1) {
798329f5518SBarry Smith       PetscReal *tmp,*tmp2;
799b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
800273d9f13SBarry Smith       tmp2 = tmp + A->N;
801273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
802064f8208SBarry Smith       *nrm = 0.0;
8033501a2bdSLois Curfman McInnes       v = mat->v;
804273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
805273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
80667e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8073501a2bdSLois Curfman McInnes         }
8083501a2bdSLois Curfman McInnes       }
809d7d1e502SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
810273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
811064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8123501a2bdSLois Curfman McInnes       }
813606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
814b0a32e0cSBarry Smith       PetscLogFlops(A->n*A->m);
8153a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
816329f5518SBarry Smith       PetscReal ntemp;
8173501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
818064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
8193a40ed3dSBarry Smith     } else {
82029bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8213501a2bdSLois Curfman McInnes     }
8223501a2bdSLois Curfman McInnes   }
8233a40ed3dSBarry Smith   PetscFunctionReturn(0);
8243501a2bdSLois Curfman McInnes }
8253501a2bdSLois Curfman McInnes 
8264a2ae208SSatish Balay #undef __FUNCT__
8274a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
8288f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8293501a2bdSLois Curfman McInnes {
8303501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8313501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8323501a2bdSLois Curfman McInnes   Mat          B;
833273d9f13SBarry Smith   int          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8343501a2bdSLois Curfman McInnes   int          j,i,ierr;
83587828ca2SBarry Smith   PetscScalar  *v;
8363501a2bdSLois Curfman McInnes 
8373a40ed3dSBarry Smith   PetscFunctionBegin;
8387c922b88SBarry Smith   if (!matout && M != N) {
83929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8407056b6fcSBarry Smith   }
841*878740d9SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B);CHKERRQ(ierr);
842*878740d9SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
843*878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
8443501a2bdSLois Curfman McInnes 
845273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
846b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr);
8473501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8483501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8493501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8503501a2bdSLois Curfman McInnes     v   += m;
8513501a2bdSLois Curfman McInnes   }
852606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8536d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8546d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8557c922b88SBarry Smith   if (matout) {
8563501a2bdSLois Curfman McInnes     *matout = B;
8573501a2bdSLois Curfman McInnes   } else {
858273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8593501a2bdSLois Curfman McInnes   }
8603a40ed3dSBarry Smith   PetscFunctionReturn(0);
861096963f5SLois Curfman McInnes }
862096963f5SLois Curfman McInnes 
863d9eff348SSatish Balay #include "petscblaslapack.h"
8644a2ae208SSatish Balay #undef __FUNCT__
8654a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
866268466fbSBarry Smith int MatScale_MPIDense(const PetscScalar *alpha,Mat inA)
86744cd7ae7SLois Curfman McInnes {
86844cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
86944cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
87044cd7ae7SLois Curfman McInnes   int          one = 1,nz;
87144cd7ae7SLois Curfman McInnes 
8723a40ed3dSBarry Smith   PetscFunctionBegin;
873273d9f13SBarry Smith   nz = inA->m*inA->N;
874268466fbSBarry Smith   BLscal_(&nz,(PetscScalar*)alpha,a->v,&one);
875b0a32e0cSBarry Smith   PetscLogFlops(nz);
8763a40ed3dSBarry Smith   PetscFunctionReturn(0);
87744cd7ae7SLois Curfman McInnes }
87844cd7ae7SLois Curfman McInnes 
8795609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
8808965ea79SLois Curfman McInnes 
8814a2ae208SSatish Balay #undef __FUNCT__
8824a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
883273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A)
884273d9f13SBarry Smith {
885273d9f13SBarry Smith   int        ierr;
886273d9f13SBarry Smith 
887273d9f13SBarry Smith   PetscFunctionBegin;
888273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
889273d9f13SBarry Smith   PetscFunctionReturn(0);
890273d9f13SBarry Smith }
891273d9f13SBarry Smith 
8928965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
89309dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
89409dc0095SBarry Smith        MatGetRow_MPIDense,
89509dc0095SBarry Smith        MatRestoreRow_MPIDense,
89609dc0095SBarry Smith        MatMult_MPIDense,
89797304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
8987c922b88SBarry Smith        MatMultTranspose_MPIDense,
8997c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9008965ea79SLois Curfman McInnes        0,
90109dc0095SBarry Smith        0,
90209dc0095SBarry Smith        0,
90397304618SKris Buschelman /*10*/ 0,
90409dc0095SBarry Smith        0,
90509dc0095SBarry Smith        0,
90609dc0095SBarry Smith        0,
90709dc0095SBarry Smith        MatTranspose_MPIDense,
90897304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
90997304618SKris Buschelman        0,
91009dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9115b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
91209dc0095SBarry Smith        MatNorm_MPIDense,
91397304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
91409dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
91509dc0095SBarry Smith        0,
91609dc0095SBarry Smith        MatSetOption_MPIDense,
91709dc0095SBarry Smith        MatZeroEntries_MPIDense,
91897304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
91909dc0095SBarry Smith        0,
92009dc0095SBarry Smith        0,
92109dc0095SBarry Smith        0,
92209dc0095SBarry Smith        0,
92397304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
924273d9f13SBarry Smith        0,
92509dc0095SBarry Smith        0,
92609dc0095SBarry Smith        MatGetArray_MPIDense,
92709dc0095SBarry Smith        MatRestoreArray_MPIDense,
92897304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
92909dc0095SBarry Smith        0,
93009dc0095SBarry Smith        0,
93109dc0095SBarry Smith        0,
93209dc0095SBarry Smith        0,
93397304618SKris Buschelman /*40*/ 0,
9342ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
93509dc0095SBarry Smith        0,
93609dc0095SBarry Smith        MatGetValues_MPIDense,
93709dc0095SBarry Smith        0,
93897304618SKris Buschelman /*45*/ 0,
93909dc0095SBarry Smith        MatScale_MPIDense,
94009dc0095SBarry Smith        0,
94109dc0095SBarry Smith        0,
94209dc0095SBarry Smith        0,
94397304618SKris Buschelman /*50*/ MatGetBlockSize_MPIDense,
94409dc0095SBarry Smith        0,
94509dc0095SBarry Smith        0,
94609dc0095SBarry Smith        0,
94709dc0095SBarry Smith        0,
94897304618SKris Buschelman /*55*/ 0,
94909dc0095SBarry Smith        0,
95009dc0095SBarry Smith        0,
95109dc0095SBarry Smith        0,
95209dc0095SBarry Smith        0,
95397304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
954b9b97703SBarry Smith        MatDestroy_MPIDense,
955b9b97703SBarry Smith        MatView_MPIDense,
95697304618SKris Buschelman        MatGetPetscMaps_Petsc,
95797304618SKris Buschelman        0,
95897304618SKris Buschelman /*65*/ 0,
95997304618SKris Buschelman        0,
96097304618SKris Buschelman        0,
96197304618SKris Buschelman        0,
96297304618SKris Buschelman        0,
96397304618SKris Buschelman /*70*/ 0,
96497304618SKris Buschelman        0,
96597304618SKris Buschelman        0,
96697304618SKris Buschelman        0,
96797304618SKris Buschelman        0,
96897304618SKris Buschelman /*75*/ 0,
96997304618SKris Buschelman        0,
97097304618SKris Buschelman        0,
97197304618SKris Buschelman        0,
97297304618SKris Buschelman        0,
97397304618SKris Buschelman /*80*/ 0,
97497304618SKris Buschelman        0,
97597304618SKris Buschelman        0,
97697304618SKris Buschelman        0,
97797304618SKris Buschelman /*85*/ MatLoad_MPIDense};
9788965ea79SLois Curfman McInnes 
979273d9f13SBarry Smith EXTERN_C_BEGIN
9804a2ae208SSatish Balay #undef __FUNCT__
981a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
982a23d5eceSKris Buschelman int MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
983a23d5eceSKris Buschelman {
984a23d5eceSKris Buschelman   Mat_MPIDense *a;
985a23d5eceSKris Buschelman   int          ierr;
986a23d5eceSKris Buschelman 
987a23d5eceSKris Buschelman   PetscFunctionBegin;
988a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
989a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
990a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
991a23d5eceSKris Buschelman 
992a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
9935c5985e7SKris Buschelman   ierr = MatCreate(PETSC_COMM_SELF,mat->m,mat->N,mat->m,mat->N,&a->A);CHKERRQ(ierr);
9945c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
9955c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
996a23d5eceSKris Buschelman   PetscLogObjectParent(mat,a->A);
997a23d5eceSKris Buschelman   PetscFunctionReturn(0);
998a23d5eceSKris Buschelman }
999a23d5eceSKris Buschelman EXTERN_C_END
1000a23d5eceSKris Buschelman 
10010bad9183SKris Buschelman /*MC
1002fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
10030bad9183SKris Buschelman 
10040bad9183SKris Buschelman    Options Database Keys:
10050bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
10060bad9183SKris Buschelman 
10070bad9183SKris Buschelman   Level: beginner
10080bad9183SKris Buschelman 
10090bad9183SKris Buschelman .seealso: MatCreateMPIDense
10100bad9183SKris Buschelman M*/
10110bad9183SKris Buschelman 
1012a23d5eceSKris Buschelman EXTERN_C_BEGIN
1013a23d5eceSKris Buschelman #undef __FUNCT__
10144a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1015273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat)
1016273d9f13SBarry Smith {
1017273d9f13SBarry Smith   Mat_MPIDense *a;
1018273d9f13SBarry Smith   int          ierr,i;
1019273d9f13SBarry Smith 
1020273d9f13SBarry Smith   PetscFunctionBegin;
1021b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1022b0a32e0cSBarry Smith   mat->data         = (void*)a;
1023273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1024273d9f13SBarry Smith   mat->factor       = 0;
1025273d9f13SBarry Smith   mat->mapping      = 0;
1026273d9f13SBarry Smith 
1027273d9f13SBarry Smith   a->factor       = 0;
1028273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1029273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
1030273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
1031273d9f13SBarry Smith 
1032273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
1033273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
1034273d9f13SBarry Smith   a->nvec = mat->n;
1035273d9f13SBarry Smith 
1036273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1037273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
10388a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
10398a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
1040273d9f13SBarry Smith 
1041273d9f13SBarry Smith   /* build local table of row and column ownerships */
1042b0a32e0cSBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1043273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
1044b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1045273d9f13SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1046273d9f13SBarry Smith   a->rowners[0] = 0;
1047273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1048273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
1049273d9f13SBarry Smith   }
1050273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
1051273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
1052273d9f13SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1053273d9f13SBarry Smith   a->cowners[0] = 0;
1054273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1055273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1056273d9f13SBarry Smith   }
1057273d9f13SBarry Smith 
1058273d9f13SBarry Smith   /* build cache for off array entries formed */
1059273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1060273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1061273d9f13SBarry Smith 
1062273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1063273d9f13SBarry Smith   a->lvec        = 0;
1064273d9f13SBarry Smith   a->Mvctx       = 0;
1065273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1066273d9f13SBarry Smith 
1067273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1068273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1069273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1070a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1071a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1072a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1073273d9f13SBarry Smith   PetscFunctionReturn(0);
1074273d9f13SBarry Smith }
1075273d9f13SBarry Smith EXTERN_C_END
1076273d9f13SBarry Smith 
1077209238afSKris Buschelman /*MC
1078002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1079209238afSKris Buschelman 
1080209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1081209238afSKris Buschelman    and MATMPIDENSE otherwise.
1082209238afSKris Buschelman 
1083209238afSKris Buschelman    Options Database Keys:
1084209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1085209238afSKris Buschelman 
1086209238afSKris Buschelman   Level: beginner
1087209238afSKris Buschelman 
1088209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1089209238afSKris Buschelman M*/
1090209238afSKris Buschelman 
1091209238afSKris Buschelman EXTERN_C_BEGIN
1092209238afSKris Buschelman #undef __FUNCT__
1093209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1094209238afSKris Buschelman int MatCreate_Dense(Mat A) {
1095209238afSKris Buschelman   int ierr,size;
1096209238afSKris Buschelman 
1097209238afSKris Buschelman   PetscFunctionBegin;
1098209238afSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr);
1099209238afSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1100209238afSKris Buschelman   if (size == 1) {
1101209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1102209238afSKris Buschelman   } else {
1103209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1104209238afSKris Buschelman   }
1105209238afSKris Buschelman   PetscFunctionReturn(0);
1106209238afSKris Buschelman }
1107209238afSKris Buschelman EXTERN_C_END
1108209238afSKris Buschelman 
11094a2ae208SSatish Balay #undef __FUNCT__
11104a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1111273d9f13SBarry Smith /*@C
1112273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1113273d9f13SBarry Smith 
1114273d9f13SBarry Smith    Not collective
1115273d9f13SBarry Smith 
1116273d9f13SBarry Smith    Input Parameters:
1117273d9f13SBarry Smith .  A - the matrix
1118273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1119273d9f13SBarry Smith    to control all matrix memory allocation.
1120273d9f13SBarry Smith 
1121273d9f13SBarry Smith    Notes:
1122273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1123273d9f13SBarry Smith    storage by columns.
1124273d9f13SBarry Smith 
1125273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1126273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1127273d9f13SBarry Smith    set data=PETSC_NULL.
1128273d9f13SBarry Smith 
1129273d9f13SBarry Smith    Level: intermediate
1130273d9f13SBarry Smith 
1131273d9f13SBarry Smith .keywords: matrix,dense, parallel
1132273d9f13SBarry Smith 
1133273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1134273d9f13SBarry Smith @*/
113587828ca2SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1136273d9f13SBarry Smith {
1137a23d5eceSKris Buschelman   int ierr,(*f)(Mat,PetscScalar *);
1138273d9f13SBarry Smith 
1139273d9f13SBarry Smith   PetscFunctionBegin;
1140565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1141a23d5eceSKris Buschelman   if (f) {
1142a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1143a23d5eceSKris Buschelman   }
1144273d9f13SBarry Smith   PetscFunctionReturn(0);
1145273d9f13SBarry Smith }
1146273d9f13SBarry Smith 
11474a2ae208SSatish Balay #undef __FUNCT__
11484a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
11498965ea79SLois Curfman McInnes /*@C
115039ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
11518965ea79SLois Curfman McInnes 
1152db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1153db81eaa0SLois Curfman McInnes 
11548965ea79SLois Curfman McInnes    Input Parameters:
1155db81eaa0SLois Curfman McInnes +  comm - MPI communicator
11568965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1157db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
11588965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1159db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
11607f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1161dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
11628965ea79SLois Curfman McInnes 
11638965ea79SLois Curfman McInnes    Output Parameter:
1164477f1c0bSLois Curfman McInnes .  A - the matrix
11658965ea79SLois Curfman McInnes 
1166b259b22eSLois Curfman McInnes    Notes:
116739ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
116839ddd567SLois Curfman McInnes    storage by columns.
11698965ea79SLois Curfman McInnes 
117018f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
117118f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
11727f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
117318f449edSLois Curfman McInnes 
11748965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
11758965ea79SLois Curfman McInnes    (possibly both).
11768965ea79SLois Curfman McInnes 
1177027ccd11SLois Curfman McInnes    Level: intermediate
1178027ccd11SLois Curfman McInnes 
117939ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
11808965ea79SLois Curfman McInnes 
118139ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
11828965ea79SLois Curfman McInnes @*/
118387828ca2SBarry Smith int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A)
11848965ea79SLois Curfman McInnes {
1185273d9f13SBarry Smith   int ierr,size;
11868965ea79SLois Curfman McInnes 
11873a40ed3dSBarry Smith   PetscFunctionBegin;
1188273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1189273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1190273d9f13SBarry Smith   if (size > 1) {
1191273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1192273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1193273d9f13SBarry Smith   } else {
1194273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1195273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
11968c469469SLois Curfman McInnes   }
11973a40ed3dSBarry Smith   PetscFunctionReturn(0);
11988965ea79SLois Curfman McInnes }
11998965ea79SLois Curfman McInnes 
12004a2ae208SSatish Balay #undef __FUNCT__
12014a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
12025609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12038965ea79SLois Curfman McInnes {
12048965ea79SLois Curfman McInnes   Mat          mat;
12053501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
120639ddd567SLois Curfman McInnes   int          ierr;
12078965ea79SLois Curfman McInnes 
12083a40ed3dSBarry Smith   PetscFunctionBegin;
12098965ea79SLois Curfman McInnes   *newmat       = 0;
1210273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1211273d9f13SBarry Smith   ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr);
1212b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1213b0a32e0cSBarry Smith   mat->data         = (void*)a;
1214549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
12153501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1216c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1217273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
12188965ea79SLois Curfman McInnes 
12198965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
12208965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
12218965ea79SLois Curfman McInnes   a->size         = oldmat->size;
12228965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1223e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1224b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
12253782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1226b0a32e0cSBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1227b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1228549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
12298798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
12308965ea79SLois Curfman McInnes 
1231329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
12325609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
1233b0a32e0cSBarry Smith   PetscLogObjectParent(mat,a->A);
12348965ea79SLois Curfman McInnes   *newmat = mat;
12353a40ed3dSBarry Smith   PetscFunctionReturn(0);
12368965ea79SLois Curfman McInnes }
12378965ea79SLois Curfman McInnes 
1238e090d566SSatish Balay #include "petscsys.h"
12398965ea79SLois Curfman McInnes 
12404a2ae208SSatish Balay #undef __FUNCT__
12414a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
124290ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat)
124390ace30eSBarry Smith {
124440011551SBarry Smith   int          *rowners,i,size,rank,m,ierr,nz,j;
124587828ca2SBarry Smith   PetscScalar  *array,*vals,*vals_ptr;
124690ace30eSBarry Smith   MPI_Status   status;
124790ace30eSBarry Smith 
12483a40ed3dSBarry Smith   PetscFunctionBegin;
1249d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1250d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
125190ace30eSBarry Smith 
125290ace30eSBarry Smith   /* determine ownership of all rows */
125390ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
1254b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1255ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
125690ace30eSBarry Smith   rowners[0] = 0;
125790ace30eSBarry Smith   for (i=2; i<=size; i++) {
125890ace30eSBarry Smith     rowners[i] += rowners[i-1];
125990ace30eSBarry Smith   }
126090ace30eSBarry Smith 
1261*878740d9SKris Buschelman   ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr);
1262*878740d9SKris Buschelman   ierr = MatSetType(*newmat,MATMPIDENSE);CHKERRQ(ierr);
1263*878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
126490ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
126590ace30eSBarry Smith 
126690ace30eSBarry Smith   if (!rank) {
126787828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
126890ace30eSBarry Smith 
126990ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
12700752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
127190ace30eSBarry Smith 
127290ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
127390ace30eSBarry Smith     vals_ptr = vals;
127490ace30eSBarry Smith     for (i=0; i<m; i++) {
127590ace30eSBarry Smith       for (j=0; j<N; j++) {
127690ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
127790ace30eSBarry Smith       }
127890ace30eSBarry Smith     }
127990ace30eSBarry Smith 
128090ace30eSBarry Smith     /* read in other processors and ship out */
128190ace30eSBarry Smith     for (i=1; i<size; i++) {
128290ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12830752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1284ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
128590ace30eSBarry Smith     }
12863a40ed3dSBarry Smith   } else {
128790ace30eSBarry Smith     /* receive numeric values */
128887828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
128990ace30eSBarry Smith 
129090ace30eSBarry Smith     /* receive message of values*/
1291ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
129290ace30eSBarry Smith 
129390ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
129490ace30eSBarry Smith     vals_ptr = vals;
129590ace30eSBarry Smith     for (i=0; i<m; i++) {
129690ace30eSBarry Smith       for (j=0; j<N; j++) {
129790ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
129890ace30eSBarry Smith       }
129990ace30eSBarry Smith     }
130090ace30eSBarry Smith   }
1301606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1302606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
13036d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13046d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13053a40ed3dSBarry Smith   PetscFunctionReturn(0);
130690ace30eSBarry Smith }
130790ace30eSBarry Smith 
13084a2ae208SSatish Balay #undef __FUNCT__
13094a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
13108e9aea5cSBarry Smith int MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat)
13118965ea79SLois Curfman McInnes {
13128965ea79SLois Curfman McInnes   Mat          A;
131387828ca2SBarry Smith   PetscScalar  *vals,*svals;
131419bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
13158965ea79SLois Curfman McInnes   MPI_Status   status;
13168965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
13178965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols;
131819bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
13193a40ed3dSBarry Smith   int          i,nz,ierr,j,rstart,rend,fd;
13208965ea79SLois Curfman McInnes 
13213a40ed3dSBarry Smith   PetscFunctionBegin;
1322d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1323d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13248965ea79SLois Curfman McInnes   if (!rank) {
1325b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
13260752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1327552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
13288965ea79SLois Curfman McInnes   }
13298965ea79SLois Curfman McInnes 
1330ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
133190ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
133290ace30eSBarry Smith 
133390ace30eSBarry Smith   /*
133490ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
133590ace30eSBarry Smith   */
133690ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
13373a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
13383a40ed3dSBarry Smith     PetscFunctionReturn(0);
133990ace30eSBarry Smith   }
134090ace30eSBarry Smith 
13418965ea79SLois Curfman McInnes   /* determine ownership of all rows */
13428965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
1343b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1344ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
13458965ea79SLois Curfman McInnes   rowners[0] = 0;
13468965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
13478965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
13488965ea79SLois Curfman McInnes   }
13498965ea79SLois Curfman McInnes   rstart = rowners[rank];
13508965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
13518965ea79SLois Curfman McInnes 
13528965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
1353b0a32e0cSBarry Smith   ierr    = PetscMalloc(2*(rend-rstart)*sizeof(int),&ourlens);CHKERRQ(ierr);
13548965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
13558965ea79SLois Curfman McInnes   if (!rank) {
1356b0a32e0cSBarry Smith     ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr);
13570752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1358b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr);
13598965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1360ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1361606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1362ca161407SBarry Smith   } else {
1363ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
13648965ea79SLois Curfman McInnes   }
13658965ea79SLois Curfman McInnes 
13668965ea79SLois Curfman McInnes   if (!rank) {
13678965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
1368b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr);
1369549d3d68SSatish Balay     ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
13708965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13718965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
13728965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
13738965ea79SLois Curfman McInnes       }
13748965ea79SLois Curfman McInnes     }
1375606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
13768965ea79SLois Curfman McInnes 
13778965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
13788965ea79SLois Curfman McInnes     maxnz = 0;
13798965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13800452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
13818965ea79SLois Curfman McInnes     }
1382b0a32e0cSBarry Smith     ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr);
13838965ea79SLois Curfman McInnes 
13848965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13858965ea79SLois Curfman McInnes     nz = procsnz[0];
1386b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13870752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13888965ea79SLois Curfman McInnes 
13898965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13908965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13918965ea79SLois Curfman McInnes       nz   = procsnz[i];
13920752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1393ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13948965ea79SLois Curfman McInnes     }
1395606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13963a40ed3dSBarry Smith   } else {
13978965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13988965ea79SLois Curfman McInnes     nz = 0;
13998965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14008965ea79SLois Curfman McInnes       nz += ourlens[i];
14018965ea79SLois Curfman McInnes     }
1402b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
14038965ea79SLois Curfman McInnes 
14048965ea79SLois Curfman McInnes     /* receive message of column indices*/
1405ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1406ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
140729bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14088965ea79SLois Curfman McInnes   }
14098965ea79SLois Curfman McInnes 
14108965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1411549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
14128965ea79SLois Curfman McInnes   jj = 0;
14138965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14148965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
14158965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
14168965ea79SLois Curfman McInnes       jj++;
14178965ea79SLois Curfman McInnes     }
14188965ea79SLois Curfman McInnes   }
14198965ea79SLois Curfman McInnes 
14208965ea79SLois Curfman McInnes   /* create our matrix */
14218965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14228965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
14238965ea79SLois Curfman McInnes   }
1424*878740d9SKris Buschelman   ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr);
1425*878740d9SKris Buschelman   ierr = MatSetType(*newmat,MATMPIDENSE);CHKERRQ(ierr);
1426*878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
14278965ea79SLois Curfman McInnes   A = *newmat;
14288965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14298965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
14308965ea79SLois Curfman McInnes   }
14318965ea79SLois Curfman McInnes 
14328965ea79SLois Curfman McInnes   if (!rank) {
143387828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14348965ea79SLois Curfman McInnes 
14358965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
14368965ea79SLois Curfman McInnes     nz = procsnz[0];
14370752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
14388965ea79SLois Curfman McInnes 
14398965ea79SLois Curfman McInnes     /* insert into matrix */
14408965ea79SLois Curfman McInnes     jj      = rstart;
14418965ea79SLois Curfman McInnes     smycols = mycols;
14428965ea79SLois Curfman McInnes     svals   = vals;
14438965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14448965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14458965ea79SLois Curfman McInnes       smycols += ourlens[i];
14468965ea79SLois Curfman McInnes       svals   += ourlens[i];
14478965ea79SLois Curfman McInnes       jj++;
14488965ea79SLois Curfman McInnes     }
14498965ea79SLois Curfman McInnes 
14508965ea79SLois Curfman McInnes     /* read in other processors and ship out */
14518965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14528965ea79SLois Curfman McInnes       nz   = procsnz[i];
14530752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1454ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
14558965ea79SLois Curfman McInnes     }
1456606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
14573a40ed3dSBarry Smith   } else {
14588965ea79SLois Curfman McInnes     /* receive numeric values */
145987828ca2SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14608965ea79SLois Curfman McInnes 
14618965ea79SLois Curfman McInnes     /* receive message of values*/
1462ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1463ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
146429bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14658965ea79SLois Curfman McInnes 
14668965ea79SLois Curfman McInnes     /* insert into matrix */
14678965ea79SLois Curfman McInnes     jj      = rstart;
14688965ea79SLois Curfman McInnes     smycols = mycols;
14698965ea79SLois Curfman McInnes     svals   = vals;
14708965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14718965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14728965ea79SLois Curfman McInnes       smycols += ourlens[i];
14738965ea79SLois Curfman McInnes       svals   += ourlens[i];
14748965ea79SLois Curfman McInnes       jj++;
14758965ea79SLois Curfman McInnes     }
14768965ea79SLois Curfman McInnes   }
1477606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1478606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1479606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1480606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
14818965ea79SLois Curfman McInnes 
14826d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14836d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14843a40ed3dSBarry Smith   PetscFunctionReturn(0);
14858965ea79SLois Curfman McInnes }
148690ace30eSBarry Smith 
148790ace30eSBarry Smith 
148890ace30eSBarry Smith 
148990ace30eSBarry Smith 
149090ace30eSBarry Smith 
1491