xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
2*be1d678aSKris Buschelman 
3ed3cc1f0SBarry Smith /*
4ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
5ed3cc1f0SBarry Smith */
6ed3cc1f0SBarry Smith 
770f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
88965ea79SLois Curfman McInnes 
9ba8c8a56SBarry Smith #undef __FUNCT__
10ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
11ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
12ba8c8a56SBarry Smith {
13ba8c8a56SBarry Smith   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
14ba8c8a56SBarry Smith   PetscErrorCode ierr;
15ba8c8a56SBarry Smith   PetscInt          lrow,rstart = mat->rstart,rend = mat->rend;
16ba8c8a56SBarry Smith 
17ba8c8a56SBarry Smith   PetscFunctionBegin;
18ba8c8a56SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
19ba8c8a56SBarry Smith   lrow = row - rstart;
20ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
21ba8c8a56SBarry Smith   PetscFunctionReturn(0);
22ba8c8a56SBarry Smith }
23ba8c8a56SBarry Smith 
24ba8c8a56SBarry Smith #undef __FUNCT__
25ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
26ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
27ba8c8a56SBarry Smith {
28ba8c8a56SBarry Smith   PetscErrorCode ierr;
29ba8c8a56SBarry Smith 
30ba8c8a56SBarry Smith   PetscFunctionBegin;
31ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
32ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
33ba8c8a56SBarry Smith   PetscFunctionReturn(0);
34ba8c8a56SBarry Smith }
35ba8c8a56SBarry Smith 
360de54da6SSatish Balay EXTERN_C_BEGIN
374a2ae208SSatish Balay #undef __FUNCT__
384a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
39*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
400de54da6SSatish Balay {
410de54da6SSatish Balay   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
426849ba73SBarry Smith   PetscErrorCode ierr;
4313f74950SBarry Smith   PetscInt          m = A->m,rstart = mdn->rstart;
4487828ca2SBarry Smith   PetscScalar  *array;
450de54da6SSatish Balay   MPI_Comm     comm;
460de54da6SSatish Balay 
470de54da6SSatish Balay   PetscFunctionBegin;
48273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
490de54da6SSatish Balay 
500de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
510de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
520de54da6SSatish Balay 
530de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
540de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
555c5985e7SKris Buschelman   ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr);
565c5985e7SKris Buschelman   ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr);
575c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
580de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
590de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
600de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
610de54da6SSatish Balay 
620de54da6SSatish Balay   *iscopy = PETSC_TRUE;
630de54da6SSatish Balay   PetscFunctionReturn(0);
640de54da6SSatish Balay }
650de54da6SSatish Balay EXTERN_C_END
660de54da6SSatish Balay 
674a2ae208SSatish Balay #undef __FUNCT__
684a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
6913f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
708965ea79SLois Curfman McInnes {
7139b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
72dfbe8321SBarry Smith   PetscErrorCode ierr;
7313f74950SBarry Smith   PetscInt          i,j,rstart = A->rstart,rend = A->rend,row;
74273d9f13SBarry Smith   PetscTruth   roworiented = A->roworiented;
758965ea79SLois Curfman McInnes 
763a40ed3dSBarry Smith   PetscFunctionBegin;
778965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
785ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
79273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
808965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
818965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
8239b7565bSBarry Smith       if (roworiented) {
8339b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
843a40ed3dSBarry Smith       } else {
858965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
865ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
87273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
8839b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
8939b7565bSBarry Smith         }
908965ea79SLois Curfman McInnes       }
913a40ed3dSBarry Smith     } else {
923782ba37SSatish Balay       if (!A->donotstash) {
9339b7565bSBarry Smith         if (roworiented) {
948798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
95d36fbae8SSatish Balay         } else {
968798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
9739b7565bSBarry Smith         }
98b49de8d1SLois Curfman McInnes       }
99b49de8d1SLois Curfman McInnes     }
1003782ba37SSatish Balay   }
1013a40ed3dSBarry Smith   PetscFunctionReturn(0);
102b49de8d1SLois Curfman McInnes }
103b49de8d1SLois Curfman McInnes 
1044a2ae208SSatish Balay #undef __FUNCT__
1054a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
10613f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
107b49de8d1SLois Curfman McInnes {
108b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
109dfbe8321SBarry Smith   PetscErrorCode ierr;
11013f74950SBarry Smith   PetscInt          i,j,rstart = mdn->rstart,rend = mdn->rend,row;
111b49de8d1SLois Curfman McInnes 
1123a40ed3dSBarry Smith   PetscFunctionBegin;
113b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
11429bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
115273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
116b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
117b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
118b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
11929bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
120273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
12129bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
122a8c6a408SBarry Smith         }
123b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
124b49de8d1SLois Curfman McInnes       }
125a8c6a408SBarry Smith     } else {
12629bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
1278965ea79SLois Curfman McInnes     }
1288965ea79SLois Curfman McInnes   }
1293a40ed3dSBarry Smith   PetscFunctionReturn(0);
1308965ea79SLois Curfman McInnes }
1318965ea79SLois Curfman McInnes 
1324a2ae208SSatish Balay #undef __FUNCT__
1334a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
134dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
135ff14e315SSatish Balay {
136ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
137dfbe8321SBarry Smith   PetscErrorCode ierr;
138ff14e315SSatish Balay 
1393a40ed3dSBarry Smith   PetscFunctionBegin;
140ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1413a40ed3dSBarry Smith   PetscFunctionReturn(0);
142ff14e315SSatish Balay }
143ff14e315SSatish Balay 
1444a2ae208SSatish Balay #undef __FUNCT__
1454a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
14613f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
147ca3fa75bSLois Curfman McInnes {
148ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
149ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
1506849ba73SBarry Smith   PetscErrorCode ierr;
15113f74950SBarry Smith   PetscInt          i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
15287828ca2SBarry Smith   PetscScalar  *av,*bv,*v = lmat->v;
153ca3fa75bSLois Curfman McInnes   Mat          newmat;
154ca3fa75bSLois Curfman McInnes 
155ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
156ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
157ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
158b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
159b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
160ca3fa75bSLois Curfman McInnes 
161ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1627eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1637eba5e9cSLois Curfman McInnes      original matrix! */
164ca3fa75bSLois Curfman McInnes 
165ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1667eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
167ca3fa75bSLois Curfman McInnes 
168ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
169ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
17029bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
1717eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
172ca3fa75bSLois Curfman McInnes     newmat = *B;
173ca3fa75bSLois Curfman McInnes   } else {
174ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
175878740d9SKris Buschelman     ierr = MatCreate(A->comm,nrows,cs,PETSC_DECIDE,ncols,&newmat);CHKERRQ(ierr);
176878740d9SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
177878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
178ca3fa75bSLois Curfman McInnes   }
179ca3fa75bSLois Curfman McInnes 
180ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
181ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
182ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
183ca3fa75bSLois Curfman McInnes 
184ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
185ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
186ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
1877eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
188ca3fa75bSLois Curfman McInnes     }
189ca3fa75bSLois Curfman McInnes   }
190ca3fa75bSLois Curfman McInnes 
191ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
192ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
193ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
194ca3fa75bSLois Curfman McInnes 
195ca3fa75bSLois Curfman McInnes   /* Free work space */
196ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
197ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
198ca3fa75bSLois Curfman McInnes   *B = newmat;
199ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
200ca3fa75bSLois Curfman McInnes }
201ca3fa75bSLois Curfman McInnes 
2024a2ae208SSatish Balay #undef __FUNCT__
2034a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
204dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
205ff14e315SSatish Balay {
2063a40ed3dSBarry Smith   PetscFunctionBegin;
2073a40ed3dSBarry Smith   PetscFunctionReturn(0);
208ff14e315SSatish Balay }
209ff14e315SSatish Balay 
2104a2ae208SSatish Balay #undef __FUNCT__
2114a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
212dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2138965ea79SLois Curfman McInnes {
21439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
2158965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
216dfbe8321SBarry Smith   PetscErrorCode ierr;
21713f74950SBarry Smith   PetscInt          nstash,reallocs;
2188965ea79SLois Curfman McInnes   InsertMode   addv;
2198965ea79SLois Curfman McInnes 
2203a40ed3dSBarry Smith   PetscFunctionBegin;
2218965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
222ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
2237056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
22429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
2258965ea79SLois Curfman McInnes   }
226e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2278965ea79SLois Curfman McInnes 
2288798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
2298798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
23063ba0a88SBarry Smith   ierr = PetscLogInfo((mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs));CHKERRQ(ierr);
2313a40ed3dSBarry Smith   PetscFunctionReturn(0);
2328965ea79SLois Curfman McInnes }
2338965ea79SLois Curfman McInnes 
2344a2ae208SSatish Balay #undef __FUNCT__
2354a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
236dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2378965ea79SLois Curfman McInnes {
23839ddd567SLois Curfman McInnes   Mat_MPIDense    *mdn=(Mat_MPIDense*)mat->data;
2396849ba73SBarry Smith   PetscErrorCode  ierr;
24013f74950SBarry Smith   PetscInt        i,*row,*col,flg,j,rstart,ncols;
24113f74950SBarry Smith   PetscMPIInt     n;
24287828ca2SBarry Smith   PetscScalar     *val;
243e0fa3b82SLois Curfman McInnes   InsertMode      addv=mat->insertmode;
2448965ea79SLois Curfman McInnes 
2453a40ed3dSBarry Smith   PetscFunctionBegin;
2468965ea79SLois Curfman McInnes   /*  wait on receives */
2477ef1d9bdSSatish Balay   while (1) {
2488798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2497ef1d9bdSSatish Balay     if (!flg) break;
2508965ea79SLois Curfman McInnes 
2517ef1d9bdSSatish Balay     for (i=0; i<n;) {
2527ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2537ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2547ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2557ef1d9bdSSatish Balay       else       ncols = n-i;
2567ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2577ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2587ef1d9bdSSatish Balay       i = j;
2598965ea79SLois Curfman McInnes     }
2607ef1d9bdSSatish Balay   }
2618798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2628965ea79SLois Curfman McInnes 
26339ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
26439ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2658965ea79SLois Curfman McInnes 
2666d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
26739ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2688965ea79SLois Curfman McInnes   }
2693a40ed3dSBarry Smith   PetscFunctionReturn(0);
2708965ea79SLois Curfman McInnes }
2718965ea79SLois Curfman McInnes 
2724a2ae208SSatish Balay #undef __FUNCT__
2734a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
274dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
2758965ea79SLois Curfman McInnes {
276dfbe8321SBarry Smith   PetscErrorCode ierr;
27739ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense*)A->data;
2783a40ed3dSBarry Smith 
2793a40ed3dSBarry Smith   PetscFunctionBegin;
2803a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2813a40ed3dSBarry Smith   PetscFunctionReturn(0);
2828965ea79SLois Curfman McInnes }
2838965ea79SLois Curfman McInnes 
2848965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2858965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2868965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2873501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2888965ea79SLois Curfman McInnes    routine.
2898965ea79SLois Curfman McInnes */
2904a2ae208SSatish Balay #undef __FUNCT__
2914a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
292dfbe8321SBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag)
2938965ea79SLois Curfman McInnes {
29439ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
2956849ba73SBarry Smith   PetscErrorCode ierr;
29613f74950SBarry Smith   PetscInt       i,N,*rows,*owners = l->rowners;
29713f74950SBarry Smith   PetscInt       *nprocs,j,idx,nsends;
29813f74950SBarry Smith   PetscInt       nmax,*svalues,*starts,*owner,nrecvs;
29913f74950SBarry Smith   PetscInt       *rvalues,tag = A->tag,count,base,slen,*source;
30013f74950SBarry Smith   PetscInt       *lens,*lrows,*values;
30113f74950SBarry Smith   PetscMPIInt    n,imdex,rank = l->rank,size = l->size;
3028965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
3038965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
3048965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
3058965ea79SLois Curfman McInnes   IS             istmp;
30635d8aa7fSBarry Smith   PetscTruth     found;
3078965ea79SLois Curfman McInnes 
3083a40ed3dSBarry Smith   PetscFunctionBegin;
309b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
3108965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
3118965ea79SLois Curfman McInnes 
3128965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
31313f74950SBarry Smith   ierr  = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
31413f74950SBarry Smith   ierr  = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
31513f74950SBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/
3168965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3178965ea79SLois Curfman McInnes     idx = rows[i];
31835d8aa7fSBarry Smith     found = PETSC_FALSE;
3198965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3208965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
321c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3228965ea79SLois Curfman McInnes       }
3238965ea79SLois Curfman McInnes     }
32429bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3258965ea79SLois Curfman McInnes   }
326c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3278965ea79SLois Curfman McInnes 
3288965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
329c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3308965ea79SLois Curfman McInnes 
3318965ea79SLois Curfman McInnes   /* post receives:   */
33213f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
333b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3348965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
33513f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3368965ea79SLois Curfman McInnes   }
3378965ea79SLois Curfman McInnes 
3388965ea79SLois Curfman McInnes   /* do sends:
3398965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3408965ea79SLois Curfman McInnes          the ith processor
3418965ea79SLois Curfman McInnes   */
34213f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
343b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
34413f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
3458965ea79SLois Curfman McInnes   starts[0]  = 0;
346c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3478965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3488965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3498965ea79SLois Curfman McInnes   }
3508965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3518965ea79SLois Curfman McInnes 
3528965ea79SLois Curfman McInnes   starts[0] = 0;
353c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3548965ea79SLois Curfman McInnes   count = 0;
3558965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
356c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
35713f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3588965ea79SLois Curfman McInnes     }
3598965ea79SLois Curfman McInnes   }
360606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3618965ea79SLois Curfman McInnes 
3628965ea79SLois Curfman McInnes   base = owners[rank];
3638965ea79SLois Curfman McInnes 
3648965ea79SLois Curfman McInnes   /*  wait on receives */
36513f74950SBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
3668965ea79SLois Curfman McInnes   source = lens + nrecvs;
3678965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3688965ea79SLois Curfman McInnes   while (count) {
369ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3708965ea79SLois Curfman McInnes     /* unpack receives into our local space */
37113f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
3728965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3738965ea79SLois Curfman McInnes     lens[imdex]    = n;
3748965ea79SLois Curfman McInnes     slen += n;
3758965ea79SLois Curfman McInnes     count--;
3768965ea79SLois Curfman McInnes   }
377606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3788965ea79SLois Curfman McInnes 
3798965ea79SLois Curfman McInnes   /* move the data into the send scatter */
38013f74950SBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
3818965ea79SLois Curfman McInnes   count = 0;
3828965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
3838965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3848965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
3858965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3868965ea79SLois Curfman McInnes     }
3878965ea79SLois Curfman McInnes   }
388606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
389606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
390606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
391606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3928965ea79SLois Curfman McInnes 
3938965ea79SLois Curfman McInnes   /* actually zap the local rows */
394029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
39552e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,istmp);CHKERRQ(ierr);
396606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3978965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3988965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3998965ea79SLois Curfman McInnes 
4008965ea79SLois Curfman McInnes   /* wait on sends */
4018965ea79SLois Curfman McInnes   if (nsends) {
402b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
403ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
404606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4058965ea79SLois Curfman McInnes   }
406606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
407606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4088965ea79SLois Curfman McInnes 
4093a40ed3dSBarry Smith   PetscFunctionReturn(0);
4108965ea79SLois Curfman McInnes }
4118965ea79SLois Curfman McInnes 
4124a2ae208SSatish Balay #undef __FUNCT__
4134a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
414dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4158965ea79SLois Curfman McInnes {
41639ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
417dfbe8321SBarry Smith   PetscErrorCode ierr;
418c456f294SBarry Smith 
4193a40ed3dSBarry Smith   PetscFunctionBegin;
42043a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
42143a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
42244cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4233a40ed3dSBarry Smith   PetscFunctionReturn(0);
4248965ea79SLois Curfman McInnes }
4258965ea79SLois Curfman McInnes 
4264a2ae208SSatish Balay #undef __FUNCT__
4274a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
428dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4298965ea79SLois Curfman McInnes {
43039ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
431dfbe8321SBarry Smith   PetscErrorCode ierr;
432c456f294SBarry Smith 
4333a40ed3dSBarry Smith   PetscFunctionBegin;
43443a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
43543a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
43644cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4373a40ed3dSBarry Smith   PetscFunctionReturn(0);
4388965ea79SLois Curfman McInnes }
4398965ea79SLois Curfman McInnes 
4404a2ae208SSatish Balay #undef __FUNCT__
4414a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
442dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
443096963f5SLois Curfman McInnes {
444096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
445dfbe8321SBarry Smith   PetscErrorCode ierr;
44687828ca2SBarry Smith   PetscScalar  zero = 0.0;
447096963f5SLois Curfman McInnes 
4483a40ed3dSBarry Smith   PetscFunctionBegin;
4493501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
4507c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
451537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
452537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4533a40ed3dSBarry Smith   PetscFunctionReturn(0);
454096963f5SLois Curfman McInnes }
455096963f5SLois Curfman McInnes 
4564a2ae208SSatish Balay #undef __FUNCT__
4574a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
458dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
459096963f5SLois Curfman McInnes {
460096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
461dfbe8321SBarry Smith   PetscErrorCode ierr;
462096963f5SLois Curfman McInnes 
4633a40ed3dSBarry Smith   PetscFunctionBegin;
4643501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4657c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
466537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
467537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4683a40ed3dSBarry Smith   PetscFunctionReturn(0);
469096963f5SLois Curfman McInnes }
470096963f5SLois Curfman McInnes 
4714a2ae208SSatish Balay #undef __FUNCT__
4724a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
473dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
4748965ea79SLois Curfman McInnes {
47539ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
476096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
477dfbe8321SBarry Smith   PetscErrorCode ierr;
47813f74950SBarry Smith   PetscInt          len,i,n,m = A->m,radd;
47987828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
480ed3cc1f0SBarry Smith 
4813a40ed3dSBarry Smith   PetscFunctionBegin;
482273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
4831ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
484096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
485273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
486273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
4877ddc982cSLois Curfman McInnes   radd = a->rstart*m;
48844cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
489096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
490096963f5SLois Curfman McInnes   }
4911ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
4923a40ed3dSBarry Smith   PetscFunctionReturn(0);
4938965ea79SLois Curfman McInnes }
4948965ea79SLois Curfman McInnes 
4954a2ae208SSatish Balay #undef __FUNCT__
4964a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
497dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
4988965ea79SLois Curfman McInnes {
4993501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
500dfbe8321SBarry Smith   PetscErrorCode ierr;
501ed3cc1f0SBarry Smith 
5023a40ed3dSBarry Smith   PetscFunctionBegin;
50394d884c6SBarry Smith 
504aa482453SBarry Smith #if defined(PETSC_USE_LOG)
50577431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N);
5068965ea79SLois Curfman McInnes #endif
5078798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
508606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
5093501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
5103501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
5113501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
512622d7880SLois Curfman McInnes   if (mdn->factor) {
513606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
514606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
515606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
516606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
517622d7880SLois Curfman McInnes   }
518606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
519901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
520901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5213a40ed3dSBarry Smith   PetscFunctionReturn(0);
5228965ea79SLois Curfman McInnes }
52339ddd567SLois Curfman McInnes 
5244a2ae208SSatish Balay #undef __FUNCT__
5254a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5266849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5278965ea79SLois Curfman McInnes {
52839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
529dfbe8321SBarry Smith   PetscErrorCode ierr;
5307056b6fcSBarry Smith 
5313a40ed3dSBarry Smith   PetscFunctionBegin;
53239ddd567SLois Curfman McInnes   if (mdn->size == 1) {
53339ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5348965ea79SLois Curfman McInnes   }
53529bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5363a40ed3dSBarry Smith   PetscFunctionReturn(0);
5378965ea79SLois Curfman McInnes }
5388965ea79SLois Curfman McInnes 
5394a2ae208SSatish Balay #undef __FUNCT__
5404a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
5416849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5428965ea79SLois Curfman McInnes {
54339ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
544dfbe8321SBarry Smith   PetscErrorCode    ierr;
54532dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
546b0a32e0cSBarry Smith   PetscViewerType   vtype;
54732077d6dSBarry Smith   PetscTruth        iascii,isdraw;
548b0a32e0cSBarry Smith   PetscViewer       sviewer;
549f3ef73ceSBarry Smith   PetscViewerFormat format;
5508965ea79SLois Curfman McInnes 
5513a40ed3dSBarry Smith   PetscFunctionBegin;
55232077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
553fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
55432077d6dSBarry Smith   if (iascii) {
555b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
556b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
557456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
5584e220ebcSLois Curfman McInnes       MatInfo info;
559888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
56077431f27SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m,
56177431f27SBarry Smith                    (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr);
562b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5633501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5643a40ed3dSBarry Smith       PetscFunctionReturn(0);
565fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5663a40ed3dSBarry Smith       PetscFunctionReturn(0);
5678965ea79SLois Curfman McInnes     }
568f1af5d2fSBarry Smith   } else if (isdraw) {
569b0a32e0cSBarry Smith     PetscDraw       draw;
570f1af5d2fSBarry Smith     PetscTruth isnull;
571f1af5d2fSBarry Smith 
572b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
573b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
574f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
575f1af5d2fSBarry Smith   }
57677ed5343SBarry Smith 
5778965ea79SLois Curfman McInnes   if (size == 1) {
57839ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5793a40ed3dSBarry Smith   } else {
5808965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5818965ea79SLois Curfman McInnes     Mat         A;
58213f74950SBarry Smith     PetscInt    M = mat->M,N = mat->N,m,row,i,nz;
583ba8c8a56SBarry Smith     PetscInt    *cols;
584ba8c8a56SBarry Smith     PetscScalar *vals;
5858965ea79SLois Curfman McInnes 
5868965ea79SLois Curfman McInnes     if (!rank) {
587878740d9SKris Buschelman       ierr = MatCreate(mat->comm,M,N,M,N,&A);CHKERRQ(ierr);
5883a40ed3dSBarry Smith     } else {
589878740d9SKris Buschelman       ierr = MatCreate(mat->comm,0,0,M,N,&A);CHKERRQ(ierr);
5908965ea79SLois Curfman McInnes     }
591be5d1d56SKris Buschelman     /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */
592878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
593878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
59452e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
5958965ea79SLois Curfman McInnes 
59639ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
59739ddd567SLois Curfman McInnes        but it's quick for now */
59851022da4SBarry Smith     A->insertmode = INSERT_VALUES;
599273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
6008965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
601ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
602ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
603ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
60439ddd567SLois Curfman McInnes       row++;
6058965ea79SLois Curfman McInnes     }
6068965ea79SLois Curfman McInnes 
6076d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6086d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
609b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
610b9b97703SBarry Smith     if (!rank) {
6116831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
6128965ea79SLois Curfman McInnes     }
613b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
614b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6158965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
6168965ea79SLois Curfman McInnes   }
6173a40ed3dSBarry Smith   PetscFunctionReturn(0);
6188965ea79SLois Curfman McInnes }
6198965ea79SLois Curfman McInnes 
6204a2ae208SSatish Balay #undef __FUNCT__
6214a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
622dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
6238965ea79SLois Curfman McInnes {
624dfbe8321SBarry Smith   PetscErrorCode ierr;
62532077d6dSBarry Smith   PetscTruth iascii,isbinary,isdraw,issocket;
6268965ea79SLois Curfman McInnes 
627433994e6SBarry Smith   PetscFunctionBegin;
6280f5bd95cSBarry Smith 
62932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
630fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
631b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
632fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6330f5bd95cSBarry Smith 
63432077d6dSBarry Smith   if (iascii || issocket || isdraw) {
635f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6360f5bd95cSBarry Smith   } else if (isbinary) {
6373a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6385cd90555SBarry Smith   } else {
639958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6408965ea79SLois Curfman McInnes   }
6413a40ed3dSBarry Smith   PetscFunctionReturn(0);
6428965ea79SLois Curfman McInnes }
6438965ea79SLois Curfman McInnes 
6444a2ae208SSatish Balay #undef __FUNCT__
6454a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
646dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6478965ea79SLois Curfman McInnes {
6483501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6493501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
650dfbe8321SBarry Smith   PetscErrorCode ierr;
651329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6528965ea79SLois Curfman McInnes 
6533a40ed3dSBarry Smith   PetscFunctionBegin;
654273d9f13SBarry Smith   info->rows_global    = (double)A->M;
655273d9f13SBarry Smith   info->columns_global = (double)A->N;
656273d9f13SBarry Smith   info->rows_local     = (double)A->m;
657273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6584e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6594e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6604e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6614e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6628965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6634e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6644e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6654e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6664e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6674e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6688965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
669d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
6704e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6714e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6724e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6734e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6744e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6758965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
676d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
6774e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6784e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6794e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6804e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6814e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6828965ea79SLois Curfman McInnes   }
6834e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6844e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6854e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6863a40ed3dSBarry Smith   PetscFunctionReturn(0);
6878965ea79SLois Curfman McInnes }
6888965ea79SLois Curfman McInnes 
6894a2ae208SSatish Balay #undef __FUNCT__
6904a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
691dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op)
6928965ea79SLois Curfman McInnes {
69339ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
694dfbe8321SBarry Smith   PetscErrorCode ierr;
6958965ea79SLois Curfman McInnes 
6963a40ed3dSBarry Smith   PetscFunctionBegin;
69712c028f9SKris Buschelman   switch (op) {
69812c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
69912c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
70012c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
70112c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
70212c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
70312c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
704273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
70512c028f9SKris Buschelman     break;
70612c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
707273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
708273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
70912c028f9SKris Buschelman     break;
71012c028f9SKris Buschelman   case MAT_ROWS_SORTED:
71112c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
71212c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
71312c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
71463ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatSetOption_MPIDense:Option ignored\n"));CHKERRQ(ierr);
71512c028f9SKris Buschelman     break;
71612c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
717273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
718273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
71912c028f9SKris Buschelman     break;
72012c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
721273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
72212c028f9SKris Buschelman     break;
72312c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
72429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
72577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
72677e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
7279a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
7289a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
7299a4540c5SBarry Smith   case MAT_HERMITIAN:
7309a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
7319a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
7329a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
73377e54ba9SKris Buschelman     break;
73412c028f9SKris Buschelman   default:
73529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
7363a40ed3dSBarry Smith   }
7373a40ed3dSBarry Smith   PetscFunctionReturn(0);
7388965ea79SLois Curfman McInnes }
7398965ea79SLois Curfman McInnes 
7408965ea79SLois Curfman McInnes 
7414a2ae208SSatish Balay #undef __FUNCT__
7424a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
743dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7445b2fa520SLois Curfman McInnes {
7455b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7465b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
74787828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
748dfbe8321SBarry Smith   PetscErrorCode ierr;
74913f74950SBarry Smith   PetscInt          i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7505b2fa520SLois Curfman McInnes 
7515b2fa520SLois Curfman McInnes   PetscFunctionBegin;
75272d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7535b2fa520SLois Curfman McInnes   if (ll) {
75472d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
75529bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
7561ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7575b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7585b2fa520SLois Curfman McInnes       x = l[i];
7595b2fa520SLois Curfman McInnes       v = mat->v + i;
7605b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7615b2fa520SLois Curfman McInnes     }
7621ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
763efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
7645b2fa520SLois Curfman McInnes   }
7655b2fa520SLois Curfman McInnes   if (rr) {
76672d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
76729bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7685b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7695b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7701ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
7715b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7725b2fa520SLois Curfman McInnes       x = r[i];
7735b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7745b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7755b2fa520SLois Curfman McInnes     }
7761ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
777efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
7785b2fa520SLois Curfman McInnes   }
7795b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7805b2fa520SLois Curfman McInnes }
7815b2fa520SLois Curfman McInnes 
7824a2ae208SSatish Balay #undef __FUNCT__
7834a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
784dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
785096963f5SLois Curfman McInnes {
7863501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7873501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
788dfbe8321SBarry Smith   PetscErrorCode ierr;
78913f74950SBarry Smith   PetscInt          i,j;
790329f5518SBarry Smith   PetscReal    sum = 0.0;
79187828ca2SBarry Smith   PetscScalar  *v = mat->v;
7923501a2bdSLois Curfman McInnes 
7933a40ed3dSBarry Smith   PetscFunctionBegin;
7943501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
795064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
7963501a2bdSLois Curfman McInnes   } else {
7973501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
798273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
799aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
800329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
8013501a2bdSLois Curfman McInnes #else
8023501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
8033501a2bdSLois Curfman McInnes #endif
8043501a2bdSLois Curfman McInnes       }
805064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
806064f8208SBarry Smith       *nrm = sqrt(*nrm);
807efee365bSSatish Balay       ierr = PetscLogFlops(2*mdn->A->n*mdn->A->m);CHKERRQ(ierr);
8083a40ed3dSBarry Smith     } else if (type == NORM_1) {
809329f5518SBarry Smith       PetscReal *tmp,*tmp2;
810b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
811273d9f13SBarry Smith       tmp2 = tmp + A->N;
812273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
813064f8208SBarry Smith       *nrm = 0.0;
8143501a2bdSLois Curfman McInnes       v = mat->v;
815273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
816273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
81767e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8183501a2bdSLois Curfman McInnes         }
8193501a2bdSLois Curfman McInnes       }
820d7d1e502SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
821273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
822064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8233501a2bdSLois Curfman McInnes       }
824606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
825efee365bSSatish Balay       ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
8263a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
827329f5518SBarry Smith       PetscReal ntemp;
8283501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
829064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
8303a40ed3dSBarry Smith     } else {
83129bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8323501a2bdSLois Curfman McInnes     }
8333501a2bdSLois Curfman McInnes   }
8343a40ed3dSBarry Smith   PetscFunctionReturn(0);
8353501a2bdSLois Curfman McInnes }
8363501a2bdSLois Curfman McInnes 
8374a2ae208SSatish Balay #undef __FUNCT__
8384a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
839dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout)
8403501a2bdSLois Curfman McInnes {
8413501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8423501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8433501a2bdSLois Curfman McInnes   Mat          B;
84413f74950SBarry Smith   PetscInt          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8456849ba73SBarry Smith   PetscErrorCode ierr;
84613f74950SBarry Smith   PetscInt          j,i;
84787828ca2SBarry Smith   PetscScalar  *v;
8483501a2bdSLois Curfman McInnes 
8493a40ed3dSBarry Smith   PetscFunctionBegin;
8507c922b88SBarry Smith   if (!matout && M != N) {
85129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8527056b6fcSBarry Smith   }
853878740d9SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B);CHKERRQ(ierr);
854878740d9SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
855878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
8563501a2bdSLois Curfman McInnes 
857273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
85813f74950SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
8593501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8603501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8613501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8623501a2bdSLois Curfman McInnes     v   += m;
8633501a2bdSLois Curfman McInnes   }
864606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8656d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8666d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8677c922b88SBarry Smith   if (matout) {
8683501a2bdSLois Curfman McInnes     *matout = B;
8693501a2bdSLois Curfman McInnes   } else {
870273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8713501a2bdSLois Curfman McInnes   }
8723a40ed3dSBarry Smith   PetscFunctionReturn(0);
873096963f5SLois Curfman McInnes }
874096963f5SLois Curfman McInnes 
875d9eff348SSatish Balay #include "petscblaslapack.h"
8764a2ae208SSatish Balay #undef __FUNCT__
8774a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
878dfbe8321SBarry Smith PetscErrorCode MatScale_MPIDense(const PetscScalar *alpha,Mat inA)
87944cd7ae7SLois Curfman McInnes {
88044cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
88144cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
8824ce68768SBarry Smith   PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N;
883efee365bSSatish Balay   PetscErrorCode ierr;
88444cd7ae7SLois Curfman McInnes 
8853a40ed3dSBarry Smith   PetscFunctionBegin;
88671044d3cSBarry Smith   BLASscal_(&nz,(PetscScalar*)alpha,a->v,&one);
887efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
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,
955521d7252SBarry Smith /*50*/ 0,
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,
989865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense,
990865e5f61SKris Buschelman        0,
991865e5f61SKris Buschelman        0,
992865e5f61SKris Buschelman        0,
993865e5f61SKris Buschelman        0,
994865e5f61SKris Buschelman        0,
995865e5f61SKris Buschelman /*90*/ 0,
996865e5f61SKris Buschelman        0,
997865e5f61SKris Buschelman        0,
998865e5f61SKris Buschelman        0,
999865e5f61SKris Buschelman        0,
1000865e5f61SKris Buschelman /*95*/ 0,
1001865e5f61SKris Buschelman        0,
1002865e5f61SKris Buschelman        0,
1003865e5f61SKris Buschelman        0};
10048965ea79SLois Curfman McInnes 
1005273d9f13SBarry Smith EXTERN_C_BEGIN
10064a2ae208SSatish Balay #undef __FUNCT__
1007a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
1008*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT 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);
102252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
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"
1041*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat)
1042273d9f13SBarry Smith {
1043273d9f13SBarry Smith   Mat_MPIDense *a;
1044dfbe8321SBarry Smith   PetscErrorCode ierr;
104513f74950SBarry Smith   PetscInt          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 */
106913f74950SBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr);
1070273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
107152e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr);
107213f74950SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_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];
107913f74950SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_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"
1121*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A)
1122dfbe8321SBarry Smith {
11236849ba73SBarry Smith   PetscErrorCode ierr;
112413f74950SBarry Smith   PetscMPIInt    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 @*/
1164*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT 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 @*/
1212*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
12138965ea79SLois Curfman McInnes {
12146849ba73SBarry Smith   PetscErrorCode ierr;
121513f74950SBarry Smith   PetscMPIInt    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;
125613f74950SBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr);
125752e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr);
125813f74950SBarry Smith   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));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);
126352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
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"
127213f74950SBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,const MatType type,Mat *newmat)
127390ace30eSBarry Smith {
12746849ba73SBarry Smith   PetscErrorCode ierr;
127532dcc486SBarry Smith   PetscMPIInt    rank,size;
127613f74950SBarry Smith   PetscInt            *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);
128613f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
128713f74950SBarry Smith   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_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;
134813f74950SBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
134913f74950SBarry Smith   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
135013f74950SBarry Smith   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
135113f74950SBarry Smith   PetscInt       i,nz,j,rstart,rend;
135213f74950SBarry Smith   int            fd;
13536849ba73SBarry Smith   PetscErrorCode ierr;
13548965ea79SLois Curfman McInnes 
13553a40ed3dSBarry Smith   PetscFunctionBegin;
1356d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1357d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13588965ea79SLois Curfman McInnes   if (!rank) {
1359b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
13600752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1361552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
13628965ea79SLois Curfman McInnes   }
13638965ea79SLois Curfman McInnes 
136413f74950SBarry Smith   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
136590ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
136690ace30eSBarry Smith 
136790ace30eSBarry Smith   /*
136890ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
136990ace30eSBarry Smith   */
137090ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1371be5d1d56SKris Buschelman     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr);
13723a40ed3dSBarry Smith     PetscFunctionReturn(0);
137390ace30eSBarry Smith   }
137490ace30eSBarry Smith 
13758965ea79SLois Curfman McInnes   /* determine ownership of all rows */
13768965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
137713f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
1378ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
13798965ea79SLois Curfman McInnes   rowners[0] = 0;
13808965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
13818965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
13828965ea79SLois Curfman McInnes   }
13838965ea79SLois Curfman McInnes   rstart = rowners[rank];
13848965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
13858965ea79SLois Curfman McInnes 
13868965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
138713f74950SBarry Smith   ierr    = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr);
13888965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
13898965ea79SLois Curfman McInnes   if (!rank) {
139013f74950SBarry Smith     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
13910752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
139213f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
13938965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
13941466f1e1SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
1395606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1396ca161407SBarry Smith   } else {
13971466f1e1SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
13988965ea79SLois Curfman McInnes   }
13998965ea79SLois Curfman McInnes 
14008965ea79SLois Curfman McInnes   if (!rank) {
14018965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
140213f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
140313f74950SBarry Smith     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
14048965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14058965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
14068965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
14078965ea79SLois Curfman McInnes       }
14088965ea79SLois Curfman McInnes     }
1409606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
14108965ea79SLois Curfman McInnes 
14118965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
14128965ea79SLois Curfman McInnes     maxnz = 0;
14138965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14140452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
14158965ea79SLois Curfman McInnes     }
141613f74950SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
14178965ea79SLois Curfman McInnes 
14188965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
14198965ea79SLois Curfman McInnes     nz = procsnz[0];
142013f74950SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
14210752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
14228965ea79SLois Curfman McInnes 
14238965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
14248965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14258965ea79SLois Curfman McInnes       nz   = procsnz[i];
14260752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
142713f74950SBarry Smith       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
14288965ea79SLois Curfman McInnes     }
1429606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
14303a40ed3dSBarry Smith   } else {
14318965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
14328965ea79SLois Curfman McInnes     nz = 0;
14338965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14348965ea79SLois Curfman McInnes       nz += ourlens[i];
14358965ea79SLois Curfman McInnes     }
143613f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
14378965ea79SLois Curfman McInnes 
14388965ea79SLois Curfman McInnes     /* receive message of column indices*/
143913f74950SBarry Smith     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
144013f74950SBarry Smith     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
144129bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14428965ea79SLois Curfman McInnes   }
14438965ea79SLois Curfman McInnes 
14448965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
144513f74950SBarry Smith   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
14468965ea79SLois Curfman McInnes   jj = 0;
14478965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14488965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
14498965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
14508965ea79SLois Curfman McInnes       jj++;
14518965ea79SLois Curfman McInnes     }
14528965ea79SLois Curfman McInnes   }
14538965ea79SLois Curfman McInnes 
14548965ea79SLois Curfman McInnes   /* create our matrix */
14558965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14568965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
14578965ea79SLois Curfman McInnes   }
1458878740d9SKris Buschelman   ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr);
1459be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1460878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
14618965ea79SLois Curfman McInnes   A = *newmat;
14628965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14638965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
14648965ea79SLois Curfman McInnes   }
14658965ea79SLois Curfman McInnes 
14668965ea79SLois Curfman McInnes   if (!rank) {
146787828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14688965ea79SLois Curfman McInnes 
14698965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
14708965ea79SLois Curfman McInnes     nz = procsnz[0];
14710752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
14728965ea79SLois Curfman McInnes 
14738965ea79SLois Curfman McInnes     /* insert into matrix */
14748965ea79SLois Curfman McInnes     jj      = rstart;
14758965ea79SLois Curfman McInnes     smycols = mycols;
14768965ea79SLois Curfman McInnes     svals   = vals;
14778965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14788965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14798965ea79SLois Curfman McInnes       smycols += ourlens[i];
14808965ea79SLois Curfman McInnes       svals   += ourlens[i];
14818965ea79SLois Curfman McInnes       jj++;
14828965ea79SLois Curfman McInnes     }
14838965ea79SLois Curfman McInnes 
14848965ea79SLois Curfman McInnes     /* read in other processors and ship out */
14858965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14868965ea79SLois Curfman McInnes       nz   = procsnz[i];
14870752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1488ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
14898965ea79SLois Curfman McInnes     }
1490606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
14913a40ed3dSBarry Smith   } else {
14928965ea79SLois Curfman McInnes     /* receive numeric values */
149384d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14948965ea79SLois Curfman McInnes 
14958965ea79SLois Curfman McInnes     /* receive message of values*/
1496ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1497ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
149829bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14998965ea79SLois Curfman McInnes 
15008965ea79SLois Curfman McInnes     /* insert into matrix */
15018965ea79SLois Curfman McInnes     jj      = rstart;
15028965ea79SLois Curfman McInnes     smycols = mycols;
15038965ea79SLois Curfman McInnes     svals   = vals;
15048965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15058965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15068965ea79SLois Curfman McInnes       smycols += ourlens[i];
15078965ea79SLois Curfman McInnes       svals   += ourlens[i];
15088965ea79SLois Curfman McInnes       jj++;
15098965ea79SLois Curfman McInnes     }
15108965ea79SLois Curfman McInnes   }
1511606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1512606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1513606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1514606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
15158965ea79SLois Curfman McInnes 
15166d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15176d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15183a40ed3dSBarry Smith   PetscFunctionReturn(0);
15198965ea79SLois Curfman McInnes }
152090ace30eSBarry Smith 
152190ace30eSBarry Smith 
152290ace30eSBarry Smith 
152390ace30eSBarry Smith 
152490ace30eSBarry Smith 
1525