xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 865e5f61f906011a1e8ca16f59ec824bec4a567e)
1ed3cc1f0SBarry Smith /*
2ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
3ed3cc1f0SBarry Smith */
4ed3cc1f0SBarry Smith 
570f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
68965ea79SLois Curfman McInnes 
70de54da6SSatish Balay EXTERN_C_BEGIN
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
10dfbe8321SBarry Smith PetscErrorCode MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
110de54da6SSatish Balay {
120de54da6SSatish Balay   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
136849ba73SBarry Smith   PetscErrorCode ierr;
146849ba73SBarry Smith   int          m = A->m,rstart = mdn->rstart;
1587828ca2SBarry Smith   PetscScalar  *array;
160de54da6SSatish Balay   MPI_Comm     comm;
170de54da6SSatish Balay 
180de54da6SSatish Balay   PetscFunctionBegin;
19273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
200de54da6SSatish Balay 
210de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
220de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
230de54da6SSatish Balay 
240de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
250de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
265c5985e7SKris Buschelman   ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr);
275c5985e7SKris Buschelman   ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr);
285c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
290de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
300de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
310de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
320de54da6SSatish Balay 
330de54da6SSatish Balay   *iscopy = PETSC_TRUE;
340de54da6SSatish Balay   PetscFunctionReturn(0);
350de54da6SSatish Balay }
360de54da6SSatish Balay EXTERN_C_END
370de54da6SSatish Balay 
384a2ae208SSatish Balay #undef __FUNCT__
394a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
40dfbe8321SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],const PetscScalar v[],InsertMode addv)
418965ea79SLois Curfman McInnes {
4239b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
43dfbe8321SBarry Smith   PetscErrorCode ierr;
44dfbe8321SBarry Smith   int          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"
77dfbe8321SBarry Smith PetscErrorCode 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;
80dfbe8321SBarry Smith   PetscErrorCode ierr;
81dfbe8321SBarry Smith   int          i,j,rstart = mdn->rstart,rend = mdn->rend,row;
82b49de8d1SLois Curfman McInnes 
833a40ed3dSBarry Smith   PetscFunctionBegin;
84b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
8529bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
86273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
87b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
88b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
89b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
9029bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
91273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
9229bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
93a8c6a408SBarry Smith         }
94b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
95b49de8d1SLois Curfman McInnes       }
96a8c6a408SBarry Smith     } else {
9729bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
988965ea79SLois Curfman McInnes     }
998965ea79SLois Curfman McInnes   }
1003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1018965ea79SLois Curfman McInnes }
1028965ea79SLois Curfman McInnes 
1034a2ae208SSatish Balay #undef __FUNCT__
1044a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
105dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
106ff14e315SSatish Balay {
107ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
108dfbe8321SBarry Smith   PetscErrorCode ierr;
109ff14e315SSatish Balay 
1103a40ed3dSBarry Smith   PetscFunctionBegin;
111ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1123a40ed3dSBarry Smith   PetscFunctionReturn(0);
113ff14e315SSatish Balay }
114ff14e315SSatish Balay 
1154a2ae208SSatish Balay #undef __FUNCT__
1164a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
1176849ba73SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
118ca3fa75bSLois Curfman McInnes {
119ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
120ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
1216849ba73SBarry Smith   PetscErrorCode ierr;
1226849ba73SBarry Smith   int          i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
12387828ca2SBarry Smith   PetscScalar  *av,*bv,*v = lmat->v;
124ca3fa75bSLois Curfman McInnes   Mat          newmat;
125ca3fa75bSLois Curfman McInnes 
126ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
127ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
128ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
129b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
130b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
131ca3fa75bSLois Curfman McInnes 
132ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1337eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1347eba5e9cSLois Curfman McInnes      original matrix! */
135ca3fa75bSLois Curfman McInnes 
136ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1377eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
138ca3fa75bSLois Curfman McInnes 
139ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
140ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
14129bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
1427eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
143ca3fa75bSLois Curfman McInnes     newmat = *B;
144ca3fa75bSLois Curfman McInnes   } else {
145ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
146878740d9SKris Buschelman     ierr = MatCreate(A->comm,nrows,cs,PETSC_DECIDE,ncols,&newmat);CHKERRQ(ierr);
147878740d9SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
148878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
149ca3fa75bSLois Curfman McInnes   }
150ca3fa75bSLois Curfman McInnes 
151ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
152ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
153ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
154ca3fa75bSLois Curfman McInnes 
155ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
156ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
157ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
1587eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
159ca3fa75bSLois Curfman McInnes     }
160ca3fa75bSLois Curfman McInnes   }
161ca3fa75bSLois Curfman McInnes 
162ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
163ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165ca3fa75bSLois Curfman McInnes 
166ca3fa75bSLois Curfman McInnes   /* Free work space */
167ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
168ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
169ca3fa75bSLois Curfman McInnes   *B = newmat;
170ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
171ca3fa75bSLois Curfman McInnes }
172ca3fa75bSLois Curfman McInnes 
1734a2ae208SSatish Balay #undef __FUNCT__
1744a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
175dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
176ff14e315SSatish Balay {
1773a40ed3dSBarry Smith   PetscFunctionBegin;
1783a40ed3dSBarry Smith   PetscFunctionReturn(0);
179ff14e315SSatish Balay }
180ff14e315SSatish Balay 
1814a2ae208SSatish Balay #undef __FUNCT__
1824a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
183dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
1848965ea79SLois Curfman McInnes {
18539ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
1868965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
187dfbe8321SBarry Smith   PetscErrorCode ierr;
188dfbe8321SBarry Smith   int          nstash,reallocs;
1898965ea79SLois Curfman McInnes   InsertMode   addv;
1908965ea79SLois Curfman McInnes 
1913a40ed3dSBarry Smith   PetscFunctionBegin;
1928965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
193ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
1947056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
19529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
1968965ea79SLois Curfman McInnes   }
197e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
1988965ea79SLois Curfman McInnes 
1998798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
2008798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
201b0a32e0cSBarry Smith   PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs);
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
2038965ea79SLois Curfman McInnes }
2048965ea79SLois Curfman McInnes 
2054a2ae208SSatish Balay #undef __FUNCT__
2064a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
207dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2088965ea79SLois Curfman McInnes {
20939ddd567SLois Curfman McInnes   Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
2106849ba73SBarry Smith   PetscErrorCode ierr;
2116849ba73SBarry Smith   int          i,n,*row,*col,flg,j,rstart,ncols;
21287828ca2SBarry Smith   PetscScalar  *val;
213e0fa3b82SLois Curfman McInnes   InsertMode   addv=mat->insertmode;
2148965ea79SLois Curfman McInnes 
2153a40ed3dSBarry Smith   PetscFunctionBegin;
2168965ea79SLois Curfman McInnes   /*  wait on receives */
2177ef1d9bdSSatish Balay   while (1) {
2188798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2197ef1d9bdSSatish Balay     if (!flg) break;
2208965ea79SLois Curfman McInnes 
2217ef1d9bdSSatish Balay     for (i=0; i<n;) {
2227ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2237ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2247ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2257ef1d9bdSSatish Balay       else       ncols = n-i;
2267ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2277ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2287ef1d9bdSSatish Balay       i = j;
2298965ea79SLois Curfman McInnes     }
2307ef1d9bdSSatish Balay   }
2318798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2328965ea79SLois Curfman McInnes 
23339ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
23439ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2358965ea79SLois Curfman McInnes 
2366d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
23739ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2388965ea79SLois Curfman McInnes   }
2393a40ed3dSBarry Smith   PetscFunctionReturn(0);
2408965ea79SLois Curfman McInnes }
2418965ea79SLois Curfman McInnes 
2424a2ae208SSatish Balay #undef __FUNCT__
2434a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
244dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
2458965ea79SLois Curfman McInnes {
246dfbe8321SBarry Smith   PetscErrorCode ierr;
24739ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense*)A->data;
2483a40ed3dSBarry Smith 
2493a40ed3dSBarry Smith   PetscFunctionBegin;
2503a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2513a40ed3dSBarry Smith   PetscFunctionReturn(0);
2528965ea79SLois Curfman McInnes }
2538965ea79SLois Curfman McInnes 
2544a2ae208SSatish Balay #undef __FUNCT__
2554a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIDense"
256dfbe8321SBarry Smith PetscErrorCode MatGetBlockSize_MPIDense(Mat A,int *bs)
2574e220ebcSLois Curfman McInnes {
2583a40ed3dSBarry Smith   PetscFunctionBegin;
2594e220ebcSLois Curfman McInnes   *bs = 1;
2603a40ed3dSBarry Smith   PetscFunctionReturn(0);
2614e220ebcSLois Curfman McInnes }
2624e220ebcSLois Curfman McInnes 
2638965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2648965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2658965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2663501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2678965ea79SLois Curfman McInnes    routine.
2688965ea79SLois Curfman McInnes */
2694a2ae208SSatish Balay #undef __FUNCT__
2704a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
271dfbe8321SBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag)
2728965ea79SLois Curfman McInnes {
27339ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
2746849ba73SBarry Smith   PetscErrorCode ierr;
2756849ba73SBarry Smith   int            i,N,*rows,*owners = l->rowners,size = l->size;
276c1dc657dSBarry Smith   int            *nprocs,j,idx,nsends;
2778965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2788965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2798965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2808965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2818965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2828965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2838965ea79SLois Curfman McInnes   IS             istmp;
28435d8aa7fSBarry Smith   PetscTruth     found;
2858965ea79SLois Curfman McInnes 
2863a40ed3dSBarry Smith   PetscFunctionBegin;
287b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
2888965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
2898965ea79SLois Curfman McInnes 
2908965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
291b0a32e0cSBarry Smith   ierr  = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr);
292549d3d68SSatish Balay   ierr  = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
293b0a32e0cSBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/
2948965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
2958965ea79SLois Curfman McInnes     idx = rows[i];
29635d8aa7fSBarry Smith     found = PETSC_FALSE;
2978965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
2988965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
299c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3008965ea79SLois Curfman McInnes       }
3018965ea79SLois Curfman McInnes     }
30229bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3038965ea79SLois Curfman McInnes   }
304c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3058965ea79SLois Curfman McInnes 
3068965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
307c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3088965ea79SLois Curfman McInnes 
3098965ea79SLois Curfman McInnes   /* post receives:   */
310b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr);
311b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3128965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
313ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3148965ea79SLois Curfman McInnes   }
3158965ea79SLois Curfman McInnes 
3168965ea79SLois Curfman McInnes   /* do sends:
3178965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3188965ea79SLois Curfman McInnes          the ith processor
3198965ea79SLois Curfman McInnes   */
320b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr);
321b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
322b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr);
3238965ea79SLois Curfman McInnes   starts[0]  = 0;
324c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3258965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3268965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3278965ea79SLois Curfman McInnes   }
3288965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3298965ea79SLois Curfman McInnes 
3308965ea79SLois Curfman McInnes   starts[0] = 0;
331c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3328965ea79SLois Curfman McInnes   count = 0;
3338965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
334c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
335c1dc657dSBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3368965ea79SLois Curfman McInnes     }
3378965ea79SLois Curfman McInnes   }
338606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3398965ea79SLois Curfman McInnes 
3408965ea79SLois Curfman McInnes   base = owners[rank];
3418965ea79SLois Curfman McInnes 
3428965ea79SLois Curfman McInnes   /*  wait on receives */
343b0a32e0cSBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr);
3448965ea79SLois Curfman McInnes   source = lens + nrecvs;
3458965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3468965ea79SLois Curfman McInnes   while (count) {
347ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3488965ea79SLois Curfman McInnes     /* unpack receives into our local space */
349ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3508965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3518965ea79SLois Curfman McInnes     lens[imdex]  = n;
3528965ea79SLois Curfman McInnes     slen += n;
3538965ea79SLois Curfman McInnes     count--;
3548965ea79SLois Curfman McInnes   }
355606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3568965ea79SLois Curfman McInnes 
3578965ea79SLois Curfman McInnes   /* move the data into the send scatter */
358b0a32e0cSBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr);
3598965ea79SLois Curfman McInnes   count = 0;
3608965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
3618965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3628965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
3638965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3648965ea79SLois Curfman McInnes     }
3658965ea79SLois Curfman McInnes   }
366606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
367606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
368606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
369606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3708965ea79SLois Curfman McInnes 
3718965ea79SLois Curfman McInnes   /* actually zap the local rows */
372029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
373b0a32e0cSBarry Smith   PetscLogObjectParent(A,istmp);
374606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3758965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3768965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes 
3788965ea79SLois Curfman McInnes   /* wait on sends */
3798965ea79SLois Curfman McInnes   if (nsends) {
380b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
381ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
382606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
3838965ea79SLois Curfman McInnes   }
384606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
385606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
3868965ea79SLois Curfman McInnes 
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
3888965ea79SLois Curfman McInnes }
3898965ea79SLois Curfman McInnes 
3904a2ae208SSatish Balay #undef __FUNCT__
3914a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
392dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3938965ea79SLois Curfman McInnes {
39439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
395dfbe8321SBarry Smith   PetscErrorCode ierr;
396c456f294SBarry Smith 
3973a40ed3dSBarry Smith   PetscFunctionBegin;
39843a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39943a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40044cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4013a40ed3dSBarry Smith   PetscFunctionReturn(0);
4028965ea79SLois Curfman McInnes }
4038965ea79SLois Curfman McInnes 
4044a2ae208SSatish Balay #undef __FUNCT__
4054a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
406dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4078965ea79SLois Curfman McInnes {
40839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
409dfbe8321SBarry Smith   PetscErrorCode ierr;
410c456f294SBarry Smith 
4113a40ed3dSBarry Smith   PetscFunctionBegin;
41243a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41343a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41444cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4153a40ed3dSBarry Smith   PetscFunctionReturn(0);
4168965ea79SLois Curfman McInnes }
4178965ea79SLois Curfman McInnes 
4184a2ae208SSatish Balay #undef __FUNCT__
4194a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
420dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
421096963f5SLois Curfman McInnes {
422096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
423dfbe8321SBarry Smith   PetscErrorCode ierr;
42487828ca2SBarry Smith   PetscScalar  zero = 0.0;
425096963f5SLois Curfman McInnes 
4263a40ed3dSBarry Smith   PetscFunctionBegin;
4273501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
4287c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
429537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
430537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4313a40ed3dSBarry Smith   PetscFunctionReturn(0);
432096963f5SLois Curfman McInnes }
433096963f5SLois Curfman McInnes 
4344a2ae208SSatish Balay #undef __FUNCT__
4354a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
436dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
437096963f5SLois Curfman McInnes {
438096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
439dfbe8321SBarry Smith   PetscErrorCode ierr;
440096963f5SLois Curfman McInnes 
4413a40ed3dSBarry Smith   PetscFunctionBegin;
4423501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4437c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
444537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
445537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4463a40ed3dSBarry Smith   PetscFunctionReturn(0);
447096963f5SLois Curfman McInnes }
448096963f5SLois Curfman McInnes 
4494a2ae208SSatish Balay #undef __FUNCT__
4504a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
451dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
4528965ea79SLois Curfman McInnes {
45339ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
454096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
455dfbe8321SBarry Smith   PetscErrorCode ierr;
456dfbe8321SBarry Smith   int          len,i,n,m = A->m,radd;
45787828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
458ed3cc1f0SBarry Smith 
4593a40ed3dSBarry Smith   PetscFunctionBegin;
460273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
4611ebc52fbSHong Zhang   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   }
4691ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
4703a40ed3dSBarry Smith   PetscFunctionReturn(0);
4718965ea79SLois Curfman McInnes }
4728965ea79SLois Curfman McInnes 
4734a2ae208SSatish Balay #undef __FUNCT__
4744a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
475dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
4768965ea79SLois Curfman McInnes {
4773501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
478dfbe8321SBarry Smith   PetscErrorCode 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 
5004a2ae208SSatish Balay #undef __FUNCT__
5014a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5026849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5038965ea79SLois Curfman McInnes {
50439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
505dfbe8321SBarry Smith   PetscErrorCode 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 
5154a2ae208SSatish Balay #undef __FUNCT__
5164a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
5176849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5188965ea79SLois Curfman McInnes {
51939ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
520dfbe8321SBarry Smith   PetscErrorCode    ierr;
52132dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
522b0a32e0cSBarry Smith   PetscViewerType   vtype;
52332077d6dSBarry Smith   PetscTruth        iascii,isdraw;
524b0a32e0cSBarry Smith   PetscViewer       sviewer;
525f3ef73ceSBarry Smith   PetscViewerFormat format;
5268965ea79SLois Curfman McInnes 
5273a40ed3dSBarry Smith   PetscFunctionBegin;
52832077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
529fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
53032077d6dSBarry Smith   if (iascii) {
531b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
532b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
533456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
5344e220ebcSLois Curfman McInnes       MatInfo info;
535888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
536b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m,
5376831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
538b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5393501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5403a40ed3dSBarry Smith       PetscFunctionReturn(0);
541fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5423a40ed3dSBarry Smith       PetscFunctionReturn(0);
5438965ea79SLois Curfman McInnes     }
544f1af5d2fSBarry Smith   } else if (isdraw) {
545b0a32e0cSBarry Smith     PetscDraw       draw;
546f1af5d2fSBarry Smith     PetscTruth isnull;
547f1af5d2fSBarry Smith 
548b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
549b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
550f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
551f1af5d2fSBarry Smith   }
55277ed5343SBarry Smith 
5538965ea79SLois Curfman McInnes   if (size == 1) {
55439ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5553a40ed3dSBarry Smith   } else {
5568965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5578965ea79SLois Curfman McInnes     Mat               A;
558b3cc6726SBarry Smith     int               M = mat->M,N = mat->N,m,row,i,nz;
559b3cc6726SBarry Smith     const int         *cols;
560b3cc6726SBarry Smith     const PetscScalar *vals;
5618965ea79SLois Curfman McInnes 
5628965ea79SLois Curfman McInnes     if (!rank) {
563878740d9SKris Buschelman       ierr = MatCreate(mat->comm,M,N,M,N,&A);CHKERRQ(ierr);
5643a40ed3dSBarry Smith     } else {
565878740d9SKris Buschelman       ierr = MatCreate(mat->comm,0,0,M,N,&A);CHKERRQ(ierr);
5668965ea79SLois Curfman McInnes     }
567be5d1d56SKris Buschelman     /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */
568878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
569878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
570b0a32e0cSBarry Smith     PetscLogObjectParent(mat,A);
5718965ea79SLois Curfman McInnes 
57239ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
57339ddd567SLois Curfman McInnes        but it's quick for now */
574273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
5758965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
57639ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57739ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
57839ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
57939ddd567SLois Curfman McInnes       row++;
5808965ea79SLois Curfman McInnes     }
5818965ea79SLois Curfman McInnes 
5826d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5836d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
584b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
585b9b97703SBarry Smith     if (!rank) {
5866831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5878965ea79SLois Curfman McInnes     }
588b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
589b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5908965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
5918965ea79SLois Curfman McInnes   }
5923a40ed3dSBarry Smith   PetscFunctionReturn(0);
5938965ea79SLois Curfman McInnes }
5948965ea79SLois Curfman McInnes 
5954a2ae208SSatish Balay #undef __FUNCT__
5964a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
597dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
5988965ea79SLois Curfman McInnes {
599dfbe8321SBarry Smith   PetscErrorCode ierr;
60032077d6dSBarry Smith   PetscTruth iascii,isbinary,isdraw,issocket;
6018965ea79SLois Curfman McInnes 
602433994e6SBarry Smith   PetscFunctionBegin;
6030f5bd95cSBarry Smith 
60432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
605fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
606b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
607fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6080f5bd95cSBarry Smith 
60932077d6dSBarry Smith   if (iascii || issocket || isdraw) {
610f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6110f5bd95cSBarry Smith   } else if (isbinary) {
6123a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6135cd90555SBarry Smith   } else {
614958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6158965ea79SLois Curfman McInnes   }
6163a40ed3dSBarry Smith   PetscFunctionReturn(0);
6178965ea79SLois Curfman McInnes }
6188965ea79SLois Curfman McInnes 
6194a2ae208SSatish Balay #undef __FUNCT__
6204a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
621dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6228965ea79SLois Curfman McInnes {
6233501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6243501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
625dfbe8321SBarry Smith   PetscErrorCode ierr;
626329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6278965ea79SLois Curfman McInnes 
6283a40ed3dSBarry Smith   PetscFunctionBegin;
629273d9f13SBarry Smith   info->rows_global    = (double)A->M;
630273d9f13SBarry Smith   info->columns_global = (double)A->N;
631273d9f13SBarry Smith   info->rows_local     = (double)A->m;
632273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6334e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6344e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6354e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6364e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6378965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6384e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6394e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6404e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6414e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6424e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6438965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
644d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
6454e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6464e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6474e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6484e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6494e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6508965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
651d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
6524e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6534e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6544e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6554e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6564e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6578965ea79SLois Curfman McInnes   }
6584e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6594e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6604e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6613a40ed3dSBarry Smith   PetscFunctionReturn(0);
6628965ea79SLois Curfman McInnes }
6638965ea79SLois Curfman McInnes 
6644a2ae208SSatish Balay #undef __FUNCT__
6654a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
666dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op)
6678965ea79SLois Curfman McInnes {
66839ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
669dfbe8321SBarry Smith   PetscErrorCode ierr;
6708965ea79SLois Curfman McInnes 
6713a40ed3dSBarry Smith   PetscFunctionBegin;
67212c028f9SKris Buschelman   switch (op) {
67312c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
67412c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
67512c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
67612c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
67712c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
67812c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
679273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
68012c028f9SKris Buschelman     break;
68112c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
682273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
683273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
68412c028f9SKris Buschelman     break;
68512c028f9SKris Buschelman   case MAT_ROWS_SORTED:
68612c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
68712c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
68812c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
689b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
69012c028f9SKris Buschelman     break;
69112c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
692273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
693273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
69412c028f9SKris Buschelman     break;
69512c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
696273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
69712c028f9SKris Buschelman     break;
69812c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
69929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
70077e54ba9SKris Buschelman   case MAT_SYMMETRIC:
70177e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
7029a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
7039a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
7049a4540c5SBarry Smith   case MAT_HERMITIAN:
7059a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
7069a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
7079a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
70877e54ba9SKris Buschelman     break;
70912c028f9SKris Buschelman   default:
71029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
7113a40ed3dSBarry Smith   }
7123a40ed3dSBarry Smith   PetscFunctionReturn(0);
7138965ea79SLois Curfman McInnes }
7148965ea79SLois Curfman McInnes 
7154a2ae208SSatish Balay #undef __FUNCT__
7164a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIDense"
717dfbe8321SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v)
7188965ea79SLois Curfman McInnes {
7193501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7206849ba73SBarry Smith   PetscErrorCode ierr;
7216849ba73SBarry Smith   int          lrow,rstart = mat->rstart,rend = mat->rend;
7228965ea79SLois Curfman McInnes 
7233a40ed3dSBarry Smith   PetscFunctionBegin;
72429bbc08cSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
7258965ea79SLois Curfman McInnes   lrow = row - rstart;
726b3cc6726SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const int **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
7273a40ed3dSBarry Smith   PetscFunctionReturn(0);
7288965ea79SLois Curfman McInnes }
7298965ea79SLois Curfman McInnes 
7304a2ae208SSatish Balay #undef __FUNCT__
7314a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIDense"
732dfbe8321SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v)
7338965ea79SLois Curfman McInnes {
734dfbe8321SBarry Smith   PetscErrorCode ierr;
735606d414cSSatish Balay 
7363a40ed3dSBarry Smith   PetscFunctionBegin;
737606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
738606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7393a40ed3dSBarry Smith   PetscFunctionReturn(0);
7408965ea79SLois Curfman McInnes }
7418965ea79SLois Curfman McInnes 
7424a2ae208SSatish Balay #undef __FUNCT__
7434a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
744dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7455b2fa520SLois Curfman McInnes {
7465b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7475b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
74887828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
749dfbe8321SBarry Smith   PetscErrorCode ierr;
750dfbe8321SBarry Smith   int          i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7515b2fa520SLois Curfman McInnes 
7525b2fa520SLois Curfman McInnes   PetscFunctionBegin;
75372d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7545b2fa520SLois Curfman McInnes   if (ll) {
75572d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
75629bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
7571ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7585b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7595b2fa520SLois Curfman McInnes       x = l[i];
7605b2fa520SLois Curfman McInnes       v = mat->v + i;
7615b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7625b2fa520SLois Curfman McInnes     }
7631ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
764b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7655b2fa520SLois Curfman McInnes   }
7665b2fa520SLois Curfman McInnes   if (rr) {
76772d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
76829bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7695b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7705b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7711ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
7725b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7735b2fa520SLois Curfman McInnes       x = r[i];
7745b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7755b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7765b2fa520SLois Curfman McInnes     }
7771ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
778b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7795b2fa520SLois Curfman McInnes   }
7805b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7815b2fa520SLois Curfman McInnes }
7825b2fa520SLois Curfman McInnes 
7834a2ae208SSatish Balay #undef __FUNCT__
7844a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
785dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
786096963f5SLois Curfman McInnes {
7873501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7883501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
789dfbe8321SBarry Smith   PetscErrorCode ierr;
790dfbe8321SBarry Smith   int          i,j;
791329f5518SBarry Smith   PetscReal    sum = 0.0;
79287828ca2SBarry Smith   PetscScalar  *v = mat->v;
7933501a2bdSLois Curfman McInnes 
7943a40ed3dSBarry Smith   PetscFunctionBegin;
7953501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
796064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
7973501a2bdSLois Curfman McInnes   } else {
7983501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
799273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
800aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
801329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
8023501a2bdSLois Curfman McInnes #else
8033501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
8043501a2bdSLois Curfman McInnes #endif
8053501a2bdSLois Curfman McInnes       }
806064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
807064f8208SBarry Smith       *nrm = sqrt(*nrm);
808b0a32e0cSBarry Smith       PetscLogFlops(2*mdn->A->n*mdn->A->m);
8093a40ed3dSBarry Smith     } else if (type == NORM_1) {
810329f5518SBarry Smith       PetscReal *tmp,*tmp2;
811b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
812273d9f13SBarry Smith       tmp2 = tmp + A->N;
813273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
814064f8208SBarry Smith       *nrm = 0.0;
8153501a2bdSLois Curfman McInnes       v = mat->v;
816273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
817273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
81867e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8193501a2bdSLois Curfman McInnes         }
8203501a2bdSLois Curfman McInnes       }
821d7d1e502SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
822273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
823064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8243501a2bdSLois Curfman McInnes       }
825606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
826b0a32e0cSBarry Smith       PetscLogFlops(A->n*A->m);
8273a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
828329f5518SBarry Smith       PetscReal ntemp;
8293501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
830064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
8313a40ed3dSBarry Smith     } else {
83229bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8333501a2bdSLois Curfman McInnes     }
8343501a2bdSLois Curfman McInnes   }
8353a40ed3dSBarry Smith   PetscFunctionReturn(0);
8363501a2bdSLois Curfman McInnes }
8373501a2bdSLois Curfman McInnes 
8384a2ae208SSatish Balay #undef __FUNCT__
8394a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
840dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout)
8413501a2bdSLois Curfman McInnes {
8423501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8433501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8443501a2bdSLois Curfman McInnes   Mat          B;
845273d9f13SBarry Smith   int          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8466849ba73SBarry Smith   PetscErrorCode ierr;
8476849ba73SBarry Smith   int          j,i;
84887828ca2SBarry Smith   PetscScalar  *v;
8493501a2bdSLois Curfman McInnes 
8503a40ed3dSBarry Smith   PetscFunctionBegin;
8517c922b88SBarry Smith   if (!matout && M != N) {
85229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8537056b6fcSBarry Smith   }
854878740d9SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B);CHKERRQ(ierr);
855878740d9SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
856878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
8573501a2bdSLois Curfman McInnes 
858273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
859b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr);
8603501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8613501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8623501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8633501a2bdSLois Curfman McInnes     v   += m;
8643501a2bdSLois Curfman McInnes   }
865606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8666d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8676d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8687c922b88SBarry Smith   if (matout) {
8693501a2bdSLois Curfman McInnes     *matout = B;
8703501a2bdSLois Curfman McInnes   } else {
871273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8723501a2bdSLois Curfman McInnes   }
8733a40ed3dSBarry Smith   PetscFunctionReturn(0);
874096963f5SLois Curfman McInnes }
875096963f5SLois Curfman McInnes 
876d9eff348SSatish Balay #include "petscblaslapack.h"
8774a2ae208SSatish Balay #undef __FUNCT__
8784a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
879dfbe8321SBarry Smith PetscErrorCode MatScale_MPIDense(const PetscScalar *alpha,Mat inA)
88044cd7ae7SLois Curfman McInnes {
88144cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
88244cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
8834ce68768SBarry Smith   PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N;
88444cd7ae7SLois Curfman McInnes 
8853a40ed3dSBarry Smith   PetscFunctionBegin;
886268466fbSBarry Smith   BLscal_(&nz,(PetscScalar*)alpha,a->v,&one);
887b0a32e0cSBarry Smith   PetscLogFlops(nz);
8883a40ed3dSBarry Smith   PetscFunctionReturn(0);
88944cd7ae7SLois Curfman McInnes }
89044cd7ae7SLois Curfman McInnes 
8916849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
8928965ea79SLois Curfman McInnes 
8934a2ae208SSatish Balay #undef __FUNCT__
8944a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
895dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
896273d9f13SBarry Smith {
897dfbe8321SBarry Smith   PetscErrorCode ierr;
898273d9f13SBarry Smith 
899273d9f13SBarry Smith   PetscFunctionBegin;
900273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
901273d9f13SBarry Smith   PetscFunctionReturn(0);
902273d9f13SBarry Smith }
903273d9f13SBarry Smith 
9048965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
90509dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
90609dc0095SBarry Smith        MatGetRow_MPIDense,
90709dc0095SBarry Smith        MatRestoreRow_MPIDense,
90809dc0095SBarry Smith        MatMult_MPIDense,
90997304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
9107c922b88SBarry Smith        MatMultTranspose_MPIDense,
9117c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9128965ea79SLois Curfman McInnes        0,
91309dc0095SBarry Smith        0,
91409dc0095SBarry Smith        0,
91597304618SKris Buschelman /*10*/ 0,
91609dc0095SBarry Smith        0,
91709dc0095SBarry Smith        0,
91809dc0095SBarry Smith        0,
91909dc0095SBarry Smith        MatTranspose_MPIDense,
92097304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
92197304618SKris Buschelman        0,
92209dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9235b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
92409dc0095SBarry Smith        MatNorm_MPIDense,
92597304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
92609dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
92709dc0095SBarry Smith        0,
92809dc0095SBarry Smith        MatSetOption_MPIDense,
92909dc0095SBarry Smith        MatZeroEntries_MPIDense,
93097304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
93109dc0095SBarry Smith        0,
93209dc0095SBarry Smith        0,
93309dc0095SBarry Smith        0,
93409dc0095SBarry Smith        0,
93597304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
936273d9f13SBarry Smith        0,
93709dc0095SBarry Smith        0,
93809dc0095SBarry Smith        MatGetArray_MPIDense,
93909dc0095SBarry Smith        MatRestoreArray_MPIDense,
94097304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
94109dc0095SBarry Smith        0,
94209dc0095SBarry Smith        0,
94309dc0095SBarry Smith        0,
94409dc0095SBarry Smith        0,
94597304618SKris Buschelman /*40*/ 0,
9462ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
94709dc0095SBarry Smith        0,
94809dc0095SBarry Smith        MatGetValues_MPIDense,
94909dc0095SBarry Smith        0,
95097304618SKris Buschelman /*45*/ 0,
95109dc0095SBarry Smith        MatScale_MPIDense,
95209dc0095SBarry Smith        0,
95309dc0095SBarry Smith        0,
95409dc0095SBarry Smith        0,
95597304618SKris Buschelman /*50*/ MatGetBlockSize_MPIDense,
95609dc0095SBarry Smith        0,
95709dc0095SBarry Smith        0,
95809dc0095SBarry Smith        0,
95909dc0095SBarry Smith        0,
96097304618SKris Buschelman /*55*/ 0,
96109dc0095SBarry Smith        0,
96209dc0095SBarry Smith        0,
96309dc0095SBarry Smith        0,
96409dc0095SBarry Smith        0,
96597304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
966b9b97703SBarry Smith        MatDestroy_MPIDense,
967b9b97703SBarry Smith        MatView_MPIDense,
96897304618SKris Buschelman        MatGetPetscMaps_Petsc,
96997304618SKris Buschelman        0,
97097304618SKris Buschelman /*65*/ 0,
97197304618SKris Buschelman        0,
97297304618SKris Buschelman        0,
97397304618SKris Buschelman        0,
97497304618SKris Buschelman        0,
97597304618SKris Buschelman /*70*/ 0,
97697304618SKris Buschelman        0,
97797304618SKris Buschelman        0,
97897304618SKris Buschelman        0,
97997304618SKris Buschelman        0,
98097304618SKris Buschelman /*75*/ 0,
98197304618SKris Buschelman        0,
98297304618SKris Buschelman        0,
98397304618SKris Buschelman        0,
98497304618SKris Buschelman        0,
98597304618SKris Buschelman /*80*/ 0,
98697304618SKris Buschelman        0,
98797304618SKris Buschelman        0,
98897304618SKris Buschelman        0,
989*865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense,
990*865e5f61SKris Buschelman        0,
991*865e5f61SKris Buschelman        0,
992*865e5f61SKris Buschelman        0,
993*865e5f61SKris Buschelman        0,
994*865e5f61SKris Buschelman        0,
995*865e5f61SKris Buschelman /*90*/ 0,
996*865e5f61SKris Buschelman        0,
997*865e5f61SKris Buschelman        0,
998*865e5f61SKris Buschelman        0,
999*865e5f61SKris Buschelman        0,
1000*865e5f61SKris Buschelman /*95*/ 0,
1001*865e5f61SKris Buschelman        0,
1002*865e5f61SKris Buschelman        0,
1003*865e5f61SKris Buschelman        0};
10048965ea79SLois Curfman McInnes 
1005273d9f13SBarry Smith EXTERN_C_BEGIN
10064a2ae208SSatish Balay #undef __FUNCT__
1007a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
1008dfbe8321SBarry Smith PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1009a23d5eceSKris Buschelman {
1010a23d5eceSKris Buschelman   Mat_MPIDense *a;
1011dfbe8321SBarry Smith   PetscErrorCode ierr;
1012a23d5eceSKris Buschelman 
1013a23d5eceSKris Buschelman   PetscFunctionBegin;
1014a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1015a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1016a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1017a23d5eceSKris Buschelman 
1018a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
10195c5985e7SKris Buschelman   ierr = MatCreate(PETSC_COMM_SELF,mat->m,mat->N,mat->m,mat->N,&a->A);CHKERRQ(ierr);
10205c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
10215c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
1022a23d5eceSKris Buschelman   PetscLogObjectParent(mat,a->A);
1023a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1024a23d5eceSKris Buschelman }
1025a23d5eceSKris Buschelman EXTERN_C_END
1026a23d5eceSKris Buschelman 
10270bad9183SKris Buschelman /*MC
1028fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
10290bad9183SKris Buschelman 
10300bad9183SKris Buschelman    Options Database Keys:
10310bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
10320bad9183SKris Buschelman 
10330bad9183SKris Buschelman   Level: beginner
10340bad9183SKris Buschelman 
10350bad9183SKris Buschelman .seealso: MatCreateMPIDense
10360bad9183SKris Buschelman M*/
10370bad9183SKris Buschelman 
1038a23d5eceSKris Buschelman EXTERN_C_BEGIN
1039a23d5eceSKris Buschelman #undef __FUNCT__
10404a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1041dfbe8321SBarry Smith PetscErrorCode MatCreate_MPIDense(Mat mat)
1042273d9f13SBarry Smith {
1043273d9f13SBarry Smith   Mat_MPIDense *a;
1044dfbe8321SBarry Smith   PetscErrorCode ierr;
1045dfbe8321SBarry Smith   int          i;
1046273d9f13SBarry Smith 
1047273d9f13SBarry Smith   PetscFunctionBegin;
1048b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1049b0a32e0cSBarry Smith   mat->data         = (void*)a;
1050273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1051273d9f13SBarry Smith   mat->factor       = 0;
1052273d9f13SBarry Smith   mat->mapping      = 0;
1053273d9f13SBarry Smith 
1054273d9f13SBarry Smith   a->factor       = 0;
1055273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1056273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
1057273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
1058273d9f13SBarry Smith 
1059273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
1060273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
1061273d9f13SBarry Smith   a->nvec = mat->n;
1062273d9f13SBarry Smith 
1063273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1064273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
10658a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
10668a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
1067273d9f13SBarry Smith 
1068273d9f13SBarry Smith   /* build local table of row and column ownerships */
1069b0a32e0cSBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1070273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
1071b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1072273d9f13SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1073273d9f13SBarry Smith   a->rowners[0] = 0;
1074273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1075273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
1076273d9f13SBarry Smith   }
1077273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
1078273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
1079273d9f13SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1080273d9f13SBarry Smith   a->cowners[0] = 0;
1081273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1082273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1083273d9f13SBarry Smith   }
1084273d9f13SBarry Smith 
1085273d9f13SBarry Smith   /* build cache for off array entries formed */
1086273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1087273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1088273d9f13SBarry Smith 
1089273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1090273d9f13SBarry Smith   a->lvec        = 0;
1091273d9f13SBarry Smith   a->Mvctx       = 0;
1092273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1093273d9f13SBarry Smith 
1094273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1095273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1096273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1097a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1098a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1099a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1100273d9f13SBarry Smith   PetscFunctionReturn(0);
1101273d9f13SBarry Smith }
1102273d9f13SBarry Smith EXTERN_C_END
1103273d9f13SBarry Smith 
1104209238afSKris Buschelman /*MC
1105002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1106209238afSKris Buschelman 
1107209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1108209238afSKris Buschelman    and MATMPIDENSE otherwise.
1109209238afSKris Buschelman 
1110209238afSKris Buschelman    Options Database Keys:
1111209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1112209238afSKris Buschelman 
1113209238afSKris Buschelman   Level: beginner
1114209238afSKris Buschelman 
1115209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1116209238afSKris Buschelman M*/
1117209238afSKris Buschelman 
1118209238afSKris Buschelman EXTERN_C_BEGIN
1119209238afSKris Buschelman #undef __FUNCT__
1120209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1121dfbe8321SBarry Smith PetscErrorCode MatCreate_Dense(Mat A)
1122dfbe8321SBarry Smith {
11236849ba73SBarry Smith   PetscErrorCode ierr;
11246849ba73SBarry Smith   int size;
1125209238afSKris Buschelman 
1126209238afSKris Buschelman   PetscFunctionBegin;
1127209238afSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr);
1128209238afSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1129209238afSKris Buschelman   if (size == 1) {
1130209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1131209238afSKris Buschelman   } else {
1132209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1133209238afSKris Buschelman   }
1134209238afSKris Buschelman   PetscFunctionReturn(0);
1135209238afSKris Buschelman }
1136209238afSKris Buschelman EXTERN_C_END
1137209238afSKris Buschelman 
11384a2ae208SSatish Balay #undef __FUNCT__
11394a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1140273d9f13SBarry Smith /*@C
1141273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1142273d9f13SBarry Smith 
1143273d9f13SBarry Smith    Not collective
1144273d9f13SBarry Smith 
1145273d9f13SBarry Smith    Input Parameters:
1146273d9f13SBarry Smith .  A - the matrix
1147273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1148273d9f13SBarry Smith    to control all matrix memory allocation.
1149273d9f13SBarry Smith 
1150273d9f13SBarry Smith    Notes:
1151273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1152273d9f13SBarry Smith    storage by columns.
1153273d9f13SBarry Smith 
1154273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1155273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1156273d9f13SBarry Smith    set data=PETSC_NULL.
1157273d9f13SBarry Smith 
1158273d9f13SBarry Smith    Level: intermediate
1159273d9f13SBarry Smith 
1160273d9f13SBarry Smith .keywords: matrix,dense, parallel
1161273d9f13SBarry Smith 
1162273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1163273d9f13SBarry Smith @*/
1164dfbe8321SBarry Smith PetscErrorCode MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1165273d9f13SBarry Smith {
1166dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar *);
1167273d9f13SBarry Smith 
1168273d9f13SBarry Smith   PetscFunctionBegin;
1169565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1170a23d5eceSKris Buschelman   if (f) {
1171a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1172a23d5eceSKris Buschelman   }
1173273d9f13SBarry Smith   PetscFunctionReturn(0);
1174273d9f13SBarry Smith }
1175273d9f13SBarry Smith 
11764a2ae208SSatish Balay #undef __FUNCT__
11774a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
11788965ea79SLois Curfman McInnes /*@C
117939ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
11808965ea79SLois Curfman McInnes 
1181db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1182db81eaa0SLois Curfman McInnes 
11838965ea79SLois Curfman McInnes    Input Parameters:
1184db81eaa0SLois Curfman McInnes +  comm - MPI communicator
11858965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1186db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
11878965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1188db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
11897f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1190dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
11918965ea79SLois Curfman McInnes 
11928965ea79SLois Curfman McInnes    Output Parameter:
1193477f1c0bSLois Curfman McInnes .  A - the matrix
11948965ea79SLois Curfman McInnes 
1195b259b22eSLois Curfman McInnes    Notes:
119639ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
119739ddd567SLois Curfman McInnes    storage by columns.
11988965ea79SLois Curfman McInnes 
119918f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
120018f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
12017f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
120218f449edSLois Curfman McInnes 
12038965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
12048965ea79SLois Curfman McInnes    (possibly both).
12058965ea79SLois Curfman McInnes 
1206027ccd11SLois Curfman McInnes    Level: intermediate
1207027ccd11SLois Curfman McInnes 
120839ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
12098965ea79SLois Curfman McInnes 
121039ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
12118965ea79SLois Curfman McInnes @*/
1212dfbe8321SBarry Smith PetscErrorCode MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A)
12138965ea79SLois Curfman McInnes {
12146849ba73SBarry Smith   PetscErrorCode ierr;
12156849ba73SBarry Smith   int size;
12168965ea79SLois Curfman McInnes 
12173a40ed3dSBarry Smith   PetscFunctionBegin;
1218273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1219273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1220273d9f13SBarry Smith   if (size > 1) {
1221273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1222273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1223273d9f13SBarry Smith   } else {
1224273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1225273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
12268c469469SLois Curfman McInnes   }
12273a40ed3dSBarry Smith   PetscFunctionReturn(0);
12288965ea79SLois Curfman McInnes }
12298965ea79SLois Curfman McInnes 
12304a2ae208SSatish Balay #undef __FUNCT__
12314a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
12326849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12338965ea79SLois Curfman McInnes {
12348965ea79SLois Curfman McInnes   Mat          mat;
12353501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
1236dfbe8321SBarry Smith   PetscErrorCode ierr;
12378965ea79SLois Curfman McInnes 
12383a40ed3dSBarry Smith   PetscFunctionBegin;
12398965ea79SLois Curfman McInnes   *newmat       = 0;
1240273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1241be5d1d56SKris Buschelman   ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr);
1242b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1243b0a32e0cSBarry Smith   mat->data         = (void*)a;
1244549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
12453501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1246c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1247273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
12488965ea79SLois Curfman McInnes 
12498965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
12508965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
12518965ea79SLois Curfman McInnes   a->size         = oldmat->size;
12528965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1253e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1254b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
12553782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1256b0a32e0cSBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1257b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1258549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
12598798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
12608965ea79SLois Curfman McInnes 
1261329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
12625609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
1263b0a32e0cSBarry Smith   PetscLogObjectParent(mat,a->A);
12648965ea79SLois Curfman McInnes   *newmat = mat;
12653a40ed3dSBarry Smith   PetscFunctionReturn(0);
12668965ea79SLois Curfman McInnes }
12678965ea79SLois Curfman McInnes 
1268e090d566SSatish Balay #include "petscsys.h"
12698965ea79SLois Curfman McInnes 
12704a2ae208SSatish Balay #undef __FUNCT__
12714a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
1272dfbe8321SBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,const MatType type,Mat *newmat)
127390ace30eSBarry Smith {
12746849ba73SBarry Smith   PetscErrorCode ierr;
127532dcc486SBarry Smith   PetscMPIInt    rank,size;
127632dcc486SBarry Smith   int            *rowners,i,m,nz,j;
127787828ca2SBarry Smith   PetscScalar    *array,*vals,*vals_ptr;
127890ace30eSBarry Smith   MPI_Status     status;
127990ace30eSBarry Smith 
12803a40ed3dSBarry Smith   PetscFunctionBegin;
1281d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1282d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
128390ace30eSBarry Smith 
128490ace30eSBarry Smith   /* determine ownership of all rows */
128590ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
1286b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1287ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
128890ace30eSBarry Smith   rowners[0] = 0;
128990ace30eSBarry Smith   for (i=2; i<=size; i++) {
129090ace30eSBarry Smith     rowners[i] += rowners[i-1];
129190ace30eSBarry Smith   }
129290ace30eSBarry Smith 
1293878740d9SKris Buschelman   ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr);
1294be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1295878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
129690ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
129790ace30eSBarry Smith 
129890ace30eSBarry Smith   if (!rank) {
129987828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
130090ace30eSBarry Smith 
130190ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
13020752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
130390ace30eSBarry Smith 
130490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
130590ace30eSBarry Smith     vals_ptr = vals;
130690ace30eSBarry Smith     for (i=0; i<m; i++) {
130790ace30eSBarry Smith       for (j=0; j<N; j++) {
130890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
130990ace30eSBarry Smith       }
131090ace30eSBarry Smith     }
131190ace30eSBarry Smith 
131290ace30eSBarry Smith     /* read in other processors and ship out */
131390ace30eSBarry Smith     for (i=1; i<size; i++) {
131490ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
13150752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1316ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
131790ace30eSBarry Smith     }
13183a40ed3dSBarry Smith   } else {
131990ace30eSBarry Smith     /* receive numeric values */
132087828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
132190ace30eSBarry Smith 
132290ace30eSBarry Smith     /* receive message of values*/
1323ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
132490ace30eSBarry Smith 
132590ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
132690ace30eSBarry Smith     vals_ptr = vals;
132790ace30eSBarry Smith     for (i=0; i<m; i++) {
132890ace30eSBarry Smith       for (j=0; j<N; j++) {
132990ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
133090ace30eSBarry Smith       }
133190ace30eSBarry Smith     }
133290ace30eSBarry Smith   }
1333606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1334606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
13356d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13366d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13373a40ed3dSBarry Smith   PetscFunctionReturn(0);
133890ace30eSBarry Smith }
133990ace30eSBarry Smith 
13404a2ae208SSatish Balay #undef __FUNCT__
13414a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
1342dfbe8321SBarry Smith PetscErrorCode MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat)
13438965ea79SLois Curfman McInnes {
13448965ea79SLois Curfman McInnes   Mat            A;
134587828ca2SBarry Smith   PetscScalar    *vals,*svals;
134619bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
13478965ea79SLois Curfman McInnes   MPI_Status     status;
134832dcc486SBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag;
134932dcc486SBarry Smith   int            header[4],*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
13508965ea79SLois Curfman McInnes   int            *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols;
13516849ba73SBarry Smith   int            i,nz,j,rstart,rend,fd;
13526849ba73SBarry Smith   PetscErrorCode ierr;
13538965ea79SLois Curfman McInnes 
13543a40ed3dSBarry Smith   PetscFunctionBegin;
1355d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1356d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13578965ea79SLois Curfman McInnes   if (!rank) {
1358b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
13590752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1360552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
13618965ea79SLois Curfman McInnes   }
13628965ea79SLois Curfman McInnes 
1363ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
136490ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
136590ace30eSBarry Smith 
136690ace30eSBarry Smith   /*
136790ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
136890ace30eSBarry Smith   */
136990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1370be5d1d56SKris Buschelman     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr);
13713a40ed3dSBarry Smith     PetscFunctionReturn(0);
137290ace30eSBarry Smith   }
137390ace30eSBarry Smith 
13748965ea79SLois Curfman McInnes   /* determine ownership of all rows */
13758965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
1376b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1377ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
13788965ea79SLois Curfman McInnes   rowners[0] = 0;
13798965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
13808965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
13818965ea79SLois Curfman McInnes   }
13828965ea79SLois Curfman McInnes   rstart = rowners[rank];
13838965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
13848965ea79SLois Curfman McInnes 
13858965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
138684d528b1SBarry Smith   ierr    = PetscMalloc(2*(rend-rstart+1)*sizeof(int),&ourlens);CHKERRQ(ierr);
13878965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
13888965ea79SLois Curfman McInnes   if (!rank) {
1389b0a32e0cSBarry Smith     ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr);
13900752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1391b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr);
13928965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1393ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1394606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1395ca161407SBarry Smith   } else {
1396ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
13978965ea79SLois Curfman McInnes   }
13988965ea79SLois Curfman McInnes 
13998965ea79SLois Curfman McInnes   if (!rank) {
14008965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
1401b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr);
1402549d3d68SSatish Balay     ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
14038965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14048965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
14058965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
14068965ea79SLois Curfman McInnes       }
14078965ea79SLois Curfman McInnes     }
1408606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
14098965ea79SLois Curfman McInnes 
14108965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
14118965ea79SLois Curfman McInnes     maxnz = 0;
14128965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14130452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
14148965ea79SLois Curfman McInnes     }
1415b0a32e0cSBarry Smith     ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr);
14168965ea79SLois Curfman McInnes 
14178965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
14188965ea79SLois Curfman McInnes     nz = procsnz[0];
1419b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
14200752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
14218965ea79SLois Curfman McInnes 
14228965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
14238965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14248965ea79SLois Curfman McInnes       nz   = procsnz[i];
14250752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1426ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
14278965ea79SLois Curfman McInnes     }
1428606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
14293a40ed3dSBarry Smith   } else {
14308965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
14318965ea79SLois Curfman McInnes     nz = 0;
14328965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14338965ea79SLois Curfman McInnes       nz += ourlens[i];
14348965ea79SLois Curfman McInnes     }
143584d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(int),&mycols);CHKERRQ(ierr);
14368965ea79SLois Curfman McInnes 
14378965ea79SLois Curfman McInnes     /* receive message of column indices*/
1438ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1439ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
144029bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14418965ea79SLois Curfman McInnes   }
14428965ea79SLois Curfman McInnes 
14438965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1444549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
14458965ea79SLois Curfman McInnes   jj = 0;
14468965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14478965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
14488965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
14498965ea79SLois Curfman McInnes       jj++;
14508965ea79SLois Curfman McInnes     }
14518965ea79SLois Curfman McInnes   }
14528965ea79SLois Curfman McInnes 
14538965ea79SLois Curfman McInnes   /* create our matrix */
14548965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14558965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
14568965ea79SLois Curfman McInnes   }
1457878740d9SKris Buschelman   ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr);
1458be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1459878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
14608965ea79SLois Curfman McInnes   A = *newmat;
14618965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14628965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
14638965ea79SLois Curfman McInnes   }
14648965ea79SLois Curfman McInnes 
14658965ea79SLois Curfman McInnes   if (!rank) {
146687828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14678965ea79SLois Curfman McInnes 
14688965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
14698965ea79SLois Curfman McInnes     nz = procsnz[0];
14700752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
14718965ea79SLois Curfman McInnes 
14728965ea79SLois Curfman McInnes     /* insert into matrix */
14738965ea79SLois Curfman McInnes     jj      = rstart;
14748965ea79SLois Curfman McInnes     smycols = mycols;
14758965ea79SLois Curfman McInnes     svals   = vals;
14768965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14778965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14788965ea79SLois Curfman McInnes       smycols += ourlens[i];
14798965ea79SLois Curfman McInnes       svals   += ourlens[i];
14808965ea79SLois Curfman McInnes       jj++;
14818965ea79SLois Curfman McInnes     }
14828965ea79SLois Curfman McInnes 
14838965ea79SLois Curfman McInnes     /* read in other processors and ship out */
14848965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14858965ea79SLois Curfman McInnes       nz   = procsnz[i];
14860752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1487ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
14888965ea79SLois Curfman McInnes     }
1489606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
14903a40ed3dSBarry Smith   } else {
14918965ea79SLois Curfman McInnes     /* receive numeric values */
149284d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14938965ea79SLois Curfman McInnes 
14948965ea79SLois Curfman McInnes     /* receive message of values*/
1495ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1496ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
149729bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14988965ea79SLois Curfman McInnes 
14998965ea79SLois Curfman McInnes     /* insert into matrix */
15008965ea79SLois Curfman McInnes     jj      = rstart;
15018965ea79SLois Curfman McInnes     smycols = mycols;
15028965ea79SLois Curfman McInnes     svals   = vals;
15038965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15048965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15058965ea79SLois Curfman McInnes       smycols += ourlens[i];
15068965ea79SLois Curfman McInnes       svals   += ourlens[i];
15078965ea79SLois Curfman McInnes       jj++;
15088965ea79SLois Curfman McInnes     }
15098965ea79SLois Curfman McInnes   }
1510606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1511606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1512606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1513606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
15148965ea79SLois Curfman McInnes 
15156d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15166d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15173a40ed3dSBarry Smith   PetscFunctionReturn(0);
15188965ea79SLois Curfman McInnes }
151990ace30eSBarry Smith 
152090ace30eSBarry Smith 
152190ace30eSBarry Smith 
152290ace30eSBarry Smith 
152390ace30eSBarry Smith 
1524