xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 52e6d16ba989af9362d0fcebb12e73dd7c9666ed)
1ed3cc1f0SBarry Smith /*
2ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
3ed3cc1f0SBarry Smith */
4ed3cc1f0SBarry Smith 
570f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
68965ea79SLois Curfman McInnes 
7ba8c8a56SBarry Smith #undef __FUNCT__
8ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
9ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
10ba8c8a56SBarry Smith {
11ba8c8a56SBarry Smith   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
12ba8c8a56SBarry Smith   PetscErrorCode ierr;
13ba8c8a56SBarry Smith   PetscInt          lrow,rstart = mat->rstart,rend = mat->rend;
14ba8c8a56SBarry Smith 
15ba8c8a56SBarry Smith   PetscFunctionBegin;
16ba8c8a56SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
17ba8c8a56SBarry Smith   lrow = row - rstart;
18ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
19ba8c8a56SBarry Smith   PetscFunctionReturn(0);
20ba8c8a56SBarry Smith }
21ba8c8a56SBarry Smith 
22ba8c8a56SBarry Smith #undef __FUNCT__
23ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
24ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
25ba8c8a56SBarry Smith {
26ba8c8a56SBarry Smith   PetscErrorCode ierr;
27ba8c8a56SBarry Smith 
28ba8c8a56SBarry Smith   PetscFunctionBegin;
29ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
30ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
31ba8c8a56SBarry Smith   PetscFunctionReturn(0);
32ba8c8a56SBarry Smith }
33ba8c8a56SBarry Smith 
340de54da6SSatish Balay EXTERN_C_BEGIN
354a2ae208SSatish Balay #undef __FUNCT__
364a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
37dfbe8321SBarry Smith PetscErrorCode MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
380de54da6SSatish Balay {
390de54da6SSatish Balay   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
406849ba73SBarry Smith   PetscErrorCode ierr;
4113f74950SBarry Smith   PetscInt          m = A->m,rstart = mdn->rstart;
4287828ca2SBarry Smith   PetscScalar  *array;
430de54da6SSatish Balay   MPI_Comm     comm;
440de54da6SSatish Balay 
450de54da6SSatish Balay   PetscFunctionBegin;
46273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
470de54da6SSatish Balay 
480de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
490de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
500de54da6SSatish Balay 
510de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
520de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
535c5985e7SKris Buschelman   ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr);
545c5985e7SKris Buschelman   ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr);
555c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
560de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
570de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
580de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
590de54da6SSatish Balay 
600de54da6SSatish Balay   *iscopy = PETSC_TRUE;
610de54da6SSatish Balay   PetscFunctionReturn(0);
620de54da6SSatish Balay }
630de54da6SSatish Balay EXTERN_C_END
640de54da6SSatish Balay 
654a2ae208SSatish Balay #undef __FUNCT__
664a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
6713f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
688965ea79SLois Curfman McInnes {
6939b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
70dfbe8321SBarry Smith   PetscErrorCode ierr;
7113f74950SBarry Smith   PetscInt          i,j,rstart = A->rstart,rend = A->rend,row;
72273d9f13SBarry Smith   PetscTruth   roworiented = A->roworiented;
738965ea79SLois Curfman McInnes 
743a40ed3dSBarry Smith   PetscFunctionBegin;
758965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
765ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
77273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
788965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
798965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
8039b7565bSBarry Smith       if (roworiented) {
8139b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
823a40ed3dSBarry Smith       } else {
838965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
845ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
85273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
8639b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
8739b7565bSBarry Smith         }
888965ea79SLois Curfman McInnes       }
893a40ed3dSBarry Smith     } else {
903782ba37SSatish Balay       if (!A->donotstash) {
9139b7565bSBarry Smith         if (roworiented) {
928798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
93d36fbae8SSatish Balay         } else {
948798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
9539b7565bSBarry Smith         }
96b49de8d1SLois Curfman McInnes       }
97b49de8d1SLois Curfman McInnes     }
983782ba37SSatish Balay   }
993a40ed3dSBarry Smith   PetscFunctionReturn(0);
100b49de8d1SLois Curfman McInnes }
101b49de8d1SLois Curfman McInnes 
1024a2ae208SSatish Balay #undef __FUNCT__
1034a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
10413f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
105b49de8d1SLois Curfman McInnes {
106b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
107dfbe8321SBarry Smith   PetscErrorCode ierr;
10813f74950SBarry Smith   PetscInt          i,j,rstart = mdn->rstart,rend = mdn->rend,row;
109b49de8d1SLois Curfman McInnes 
1103a40ed3dSBarry Smith   PetscFunctionBegin;
111b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
11229bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
113273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
114b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
115b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
116b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
11729bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
118273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
11929bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
120a8c6a408SBarry Smith         }
121b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
122b49de8d1SLois Curfman McInnes       }
123a8c6a408SBarry Smith     } else {
12429bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
1258965ea79SLois Curfman McInnes     }
1268965ea79SLois Curfman McInnes   }
1273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1288965ea79SLois Curfman McInnes }
1298965ea79SLois Curfman McInnes 
1304a2ae208SSatish Balay #undef __FUNCT__
1314a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
132dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
133ff14e315SSatish Balay {
134ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
135dfbe8321SBarry Smith   PetscErrorCode ierr;
136ff14e315SSatish Balay 
1373a40ed3dSBarry Smith   PetscFunctionBegin;
138ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1393a40ed3dSBarry Smith   PetscFunctionReturn(0);
140ff14e315SSatish Balay }
141ff14e315SSatish Balay 
1424a2ae208SSatish Balay #undef __FUNCT__
1434a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
14413f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
145ca3fa75bSLois Curfman McInnes {
146ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
147ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
1486849ba73SBarry Smith   PetscErrorCode ierr;
14913f74950SBarry Smith   PetscInt          i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
15087828ca2SBarry Smith   PetscScalar  *av,*bv,*v = lmat->v;
151ca3fa75bSLois Curfman McInnes   Mat          newmat;
152ca3fa75bSLois Curfman McInnes 
153ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
154ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
155ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
156b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
157b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
158ca3fa75bSLois Curfman McInnes 
159ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1607eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1617eba5e9cSLois Curfman McInnes      original matrix! */
162ca3fa75bSLois Curfman McInnes 
163ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1647eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
165ca3fa75bSLois Curfman McInnes 
166ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
167ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
16829bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
1697eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
170ca3fa75bSLois Curfman McInnes     newmat = *B;
171ca3fa75bSLois Curfman McInnes   } else {
172ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
173878740d9SKris Buschelman     ierr = MatCreate(A->comm,nrows,cs,PETSC_DECIDE,ncols,&newmat);CHKERRQ(ierr);
174878740d9SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
175878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
176ca3fa75bSLois Curfman McInnes   }
177ca3fa75bSLois Curfman McInnes 
178ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
179ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
180ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
181ca3fa75bSLois Curfman McInnes 
182ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
183ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
184ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
1857eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
186ca3fa75bSLois Curfman McInnes     }
187ca3fa75bSLois Curfman McInnes   }
188ca3fa75bSLois Curfman McInnes 
189ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
190ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
191ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
192ca3fa75bSLois Curfman McInnes 
193ca3fa75bSLois Curfman McInnes   /* Free work space */
194ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
195ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
196ca3fa75bSLois Curfman McInnes   *B = newmat;
197ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
198ca3fa75bSLois Curfman McInnes }
199ca3fa75bSLois Curfman McInnes 
2004a2ae208SSatish Balay #undef __FUNCT__
2014a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
202dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
203ff14e315SSatish Balay {
2043a40ed3dSBarry Smith   PetscFunctionBegin;
2053a40ed3dSBarry Smith   PetscFunctionReturn(0);
206ff14e315SSatish Balay }
207ff14e315SSatish Balay 
2084a2ae208SSatish Balay #undef __FUNCT__
2094a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
210dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2118965ea79SLois Curfman McInnes {
21239ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
2138965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
214dfbe8321SBarry Smith   PetscErrorCode ierr;
21513f74950SBarry Smith   PetscInt          nstash,reallocs;
2168965ea79SLois Curfman McInnes   InsertMode   addv;
2178965ea79SLois Curfman McInnes 
2183a40ed3dSBarry Smith   PetscFunctionBegin;
2198965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
220ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
2217056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
22229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
2238965ea79SLois Curfman McInnes   }
224e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2258965ea79SLois Curfman McInnes 
2268798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
2278798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
22877431f27SBarry Smith   PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
2293a40ed3dSBarry Smith   PetscFunctionReturn(0);
2308965ea79SLois Curfman McInnes }
2318965ea79SLois Curfman McInnes 
2324a2ae208SSatish Balay #undef __FUNCT__
2334a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
234dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2358965ea79SLois Curfman McInnes {
23639ddd567SLois Curfman McInnes   Mat_MPIDense    *mdn=(Mat_MPIDense*)mat->data;
2376849ba73SBarry Smith   PetscErrorCode  ierr;
23813f74950SBarry Smith   PetscInt        i,*row,*col,flg,j,rstart,ncols;
23913f74950SBarry Smith   PetscMPIInt     n;
24087828ca2SBarry Smith   PetscScalar     *val;
241e0fa3b82SLois Curfman McInnes   InsertMode      addv=mat->insertmode;
2428965ea79SLois Curfman McInnes 
2433a40ed3dSBarry Smith   PetscFunctionBegin;
2448965ea79SLois Curfman McInnes   /*  wait on receives */
2457ef1d9bdSSatish Balay   while (1) {
2468798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2477ef1d9bdSSatish Balay     if (!flg) break;
2488965ea79SLois Curfman McInnes 
2497ef1d9bdSSatish Balay     for (i=0; i<n;) {
2507ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2517ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2527ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2537ef1d9bdSSatish Balay       else       ncols = n-i;
2547ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2557ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2567ef1d9bdSSatish Balay       i = j;
2578965ea79SLois Curfman McInnes     }
2587ef1d9bdSSatish Balay   }
2598798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2608965ea79SLois Curfman McInnes 
26139ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
26239ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2638965ea79SLois Curfman McInnes 
2646d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
26539ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2668965ea79SLois Curfman McInnes   }
2673a40ed3dSBarry Smith   PetscFunctionReturn(0);
2688965ea79SLois Curfman McInnes }
2698965ea79SLois Curfman McInnes 
2704a2ae208SSatish Balay #undef __FUNCT__
2714a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
272dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
2738965ea79SLois Curfman McInnes {
274dfbe8321SBarry Smith   PetscErrorCode ierr;
27539ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense*)A->data;
2763a40ed3dSBarry Smith 
2773a40ed3dSBarry Smith   PetscFunctionBegin;
2783a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2793a40ed3dSBarry Smith   PetscFunctionReturn(0);
2808965ea79SLois Curfman McInnes }
2818965ea79SLois Curfman McInnes 
2828965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2838965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2848965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2853501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2868965ea79SLois Curfman McInnes    routine.
2878965ea79SLois Curfman McInnes */
2884a2ae208SSatish Balay #undef __FUNCT__
2894a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
290dfbe8321SBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag)
2918965ea79SLois Curfman McInnes {
29239ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
2936849ba73SBarry Smith   PetscErrorCode ierr;
29413f74950SBarry Smith   PetscInt       i,N,*rows,*owners = l->rowners;
29513f74950SBarry Smith   PetscInt       *nprocs,j,idx,nsends;
29613f74950SBarry Smith   PetscInt       nmax,*svalues,*starts,*owner,nrecvs;
29713f74950SBarry Smith   PetscInt       *rvalues,tag = A->tag,count,base,slen,*source;
29813f74950SBarry Smith   PetscInt       *lens,*lrows,*values;
29913f74950SBarry Smith   PetscMPIInt    n,imdex,rank = l->rank,size = l->size;
3008965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
3018965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
3028965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
3038965ea79SLois Curfman McInnes   IS             istmp;
30435d8aa7fSBarry Smith   PetscTruth     found;
3058965ea79SLois Curfman McInnes 
3063a40ed3dSBarry Smith   PetscFunctionBegin;
307b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
3088965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
3098965ea79SLois Curfman McInnes 
3108965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
31113f74950SBarry Smith   ierr  = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
31213f74950SBarry Smith   ierr  = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
31313f74950SBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/
3148965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3158965ea79SLois Curfman McInnes     idx = rows[i];
31635d8aa7fSBarry Smith     found = PETSC_FALSE;
3178965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3188965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
319c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3208965ea79SLois Curfman McInnes       }
3218965ea79SLois Curfman McInnes     }
32229bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3238965ea79SLois Curfman McInnes   }
324c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3258965ea79SLois Curfman McInnes 
3268965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
327c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3288965ea79SLois Curfman McInnes 
3298965ea79SLois Curfman McInnes   /* post receives:   */
33013f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
331b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3328965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
33313f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3348965ea79SLois Curfman McInnes   }
3358965ea79SLois Curfman McInnes 
3368965ea79SLois Curfman McInnes   /* do sends:
3378965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3388965ea79SLois Curfman McInnes          the ith processor
3398965ea79SLois Curfman McInnes   */
34013f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
341b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
34213f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
3438965ea79SLois Curfman McInnes   starts[0]  = 0;
344c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3458965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3468965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3478965ea79SLois Curfman McInnes   }
3488965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3498965ea79SLois Curfman McInnes 
3508965ea79SLois Curfman McInnes   starts[0] = 0;
351c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3528965ea79SLois Curfman McInnes   count = 0;
3538965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
354c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
35513f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3568965ea79SLois Curfman McInnes     }
3578965ea79SLois Curfman McInnes   }
358606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3598965ea79SLois Curfman McInnes 
3608965ea79SLois Curfman McInnes   base = owners[rank];
3618965ea79SLois Curfman McInnes 
3628965ea79SLois Curfman McInnes   /*  wait on receives */
36313f74950SBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
3648965ea79SLois Curfman McInnes   source = lens + nrecvs;
3658965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3668965ea79SLois Curfman McInnes   while (count) {
367ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3688965ea79SLois Curfman McInnes     /* unpack receives into our local space */
36913f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
3708965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3718965ea79SLois Curfman McInnes     lens[imdex]    = n;
3728965ea79SLois Curfman McInnes     slen += n;
3738965ea79SLois Curfman McInnes     count--;
3748965ea79SLois Curfman McInnes   }
375606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3768965ea79SLois Curfman McInnes 
3778965ea79SLois Curfman McInnes   /* move the data into the send scatter */
37813f74950SBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
3798965ea79SLois Curfman McInnes   count = 0;
3808965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
3818965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3828965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
3838965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3848965ea79SLois Curfman McInnes     }
3858965ea79SLois Curfman McInnes   }
386606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
387606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
388606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
389606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3908965ea79SLois Curfman McInnes 
3918965ea79SLois Curfman McInnes   /* actually zap the local rows */
392029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
393*52e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,istmp);CHKERRQ(ierr);
394606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3958965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3968965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3978965ea79SLois Curfman McInnes 
3988965ea79SLois Curfman McInnes   /* wait on sends */
3998965ea79SLois Curfman McInnes   if (nsends) {
400b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
401ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
402606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4038965ea79SLois Curfman McInnes   }
404606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
405606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4068965ea79SLois Curfman McInnes 
4073a40ed3dSBarry Smith   PetscFunctionReturn(0);
4088965ea79SLois Curfman McInnes }
4098965ea79SLois Curfman McInnes 
4104a2ae208SSatish Balay #undef __FUNCT__
4114a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
412dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4138965ea79SLois Curfman McInnes {
41439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
415dfbe8321SBarry Smith   PetscErrorCode ierr;
416c456f294SBarry Smith 
4173a40ed3dSBarry Smith   PetscFunctionBegin;
41843a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41943a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
42044cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4213a40ed3dSBarry Smith   PetscFunctionReturn(0);
4228965ea79SLois Curfman McInnes }
4238965ea79SLois Curfman McInnes 
4244a2ae208SSatish Balay #undef __FUNCT__
4254a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
426dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4278965ea79SLois Curfman McInnes {
42839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
429dfbe8321SBarry Smith   PetscErrorCode ierr;
430c456f294SBarry Smith 
4313a40ed3dSBarry Smith   PetscFunctionBegin;
43243a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
43343a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
43444cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4353a40ed3dSBarry Smith   PetscFunctionReturn(0);
4368965ea79SLois Curfman McInnes }
4378965ea79SLois Curfman McInnes 
4384a2ae208SSatish Balay #undef __FUNCT__
4394a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
440dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
441096963f5SLois Curfman McInnes {
442096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
443dfbe8321SBarry Smith   PetscErrorCode ierr;
44487828ca2SBarry Smith   PetscScalar  zero = 0.0;
445096963f5SLois Curfman McInnes 
4463a40ed3dSBarry Smith   PetscFunctionBegin;
4473501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
4487c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
449537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
450537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4513a40ed3dSBarry Smith   PetscFunctionReturn(0);
452096963f5SLois Curfman McInnes }
453096963f5SLois Curfman McInnes 
4544a2ae208SSatish Balay #undef __FUNCT__
4554a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
456dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
457096963f5SLois Curfman McInnes {
458096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
459dfbe8321SBarry Smith   PetscErrorCode ierr;
460096963f5SLois Curfman McInnes 
4613a40ed3dSBarry Smith   PetscFunctionBegin;
4623501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4637c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
464537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
465537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4663a40ed3dSBarry Smith   PetscFunctionReturn(0);
467096963f5SLois Curfman McInnes }
468096963f5SLois Curfman McInnes 
4694a2ae208SSatish Balay #undef __FUNCT__
4704a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
471dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
4728965ea79SLois Curfman McInnes {
47339ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
474096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
475dfbe8321SBarry Smith   PetscErrorCode ierr;
47613f74950SBarry Smith   PetscInt          len,i,n,m = A->m,radd;
47787828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
478ed3cc1f0SBarry Smith 
4793a40ed3dSBarry Smith   PetscFunctionBegin;
480273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
4811ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
482096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
483273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
484273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
4857ddc982cSLois Curfman McInnes   radd = a->rstart*m;
48644cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
487096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
488096963f5SLois Curfman McInnes   }
4891ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
4903a40ed3dSBarry Smith   PetscFunctionReturn(0);
4918965ea79SLois Curfman McInnes }
4928965ea79SLois Curfman McInnes 
4934a2ae208SSatish Balay #undef __FUNCT__
4944a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
495dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
4968965ea79SLois Curfman McInnes {
4973501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
498dfbe8321SBarry Smith   PetscErrorCode ierr;
499ed3cc1f0SBarry Smith 
5003a40ed3dSBarry Smith   PetscFunctionBegin;
50194d884c6SBarry Smith 
502aa482453SBarry Smith #if defined(PETSC_USE_LOG)
50377431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N);
5048965ea79SLois Curfman McInnes #endif
5058798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
506606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
5073501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
5083501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
5093501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
510622d7880SLois Curfman McInnes   if (mdn->factor) {
511606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
512606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
513606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
514606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
515622d7880SLois Curfman McInnes   }
516606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
517901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
518901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5193a40ed3dSBarry Smith   PetscFunctionReturn(0);
5208965ea79SLois Curfman McInnes }
52139ddd567SLois Curfman McInnes 
5224a2ae208SSatish Balay #undef __FUNCT__
5234a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5246849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5258965ea79SLois Curfman McInnes {
52639ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
527dfbe8321SBarry Smith   PetscErrorCode ierr;
5287056b6fcSBarry Smith 
5293a40ed3dSBarry Smith   PetscFunctionBegin;
53039ddd567SLois Curfman McInnes   if (mdn->size == 1) {
53139ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5328965ea79SLois Curfman McInnes   }
53329bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5343a40ed3dSBarry Smith   PetscFunctionReturn(0);
5358965ea79SLois Curfman McInnes }
5368965ea79SLois Curfman McInnes 
5374a2ae208SSatish Balay #undef __FUNCT__
5384a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
5396849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5408965ea79SLois Curfman McInnes {
54139ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
542dfbe8321SBarry Smith   PetscErrorCode    ierr;
54332dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
544b0a32e0cSBarry Smith   PetscViewerType   vtype;
54532077d6dSBarry Smith   PetscTruth        iascii,isdraw;
546b0a32e0cSBarry Smith   PetscViewer       sviewer;
547f3ef73ceSBarry Smith   PetscViewerFormat format;
5488965ea79SLois Curfman McInnes 
5493a40ed3dSBarry Smith   PetscFunctionBegin;
55032077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
551fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
55232077d6dSBarry Smith   if (iascii) {
553b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
554b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
555456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
5564e220ebcSLois Curfman McInnes       MatInfo info;
557888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
55877431f27SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m,
55977431f27SBarry Smith                    (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr);
560b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5613501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5623a40ed3dSBarry Smith       PetscFunctionReturn(0);
563fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5643a40ed3dSBarry Smith       PetscFunctionReturn(0);
5658965ea79SLois Curfman McInnes     }
566f1af5d2fSBarry Smith   } else if (isdraw) {
567b0a32e0cSBarry Smith     PetscDraw       draw;
568f1af5d2fSBarry Smith     PetscTruth isnull;
569f1af5d2fSBarry Smith 
570b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
571b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
572f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
573f1af5d2fSBarry Smith   }
57477ed5343SBarry Smith 
5758965ea79SLois Curfman McInnes   if (size == 1) {
57639ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5773a40ed3dSBarry Smith   } else {
5788965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5798965ea79SLois Curfman McInnes     Mat         A;
58013f74950SBarry Smith     PetscInt    M = mat->M,N = mat->N,m,row,i,nz;
581ba8c8a56SBarry Smith     PetscInt    *cols;
582ba8c8a56SBarry Smith     PetscScalar *vals;
5838965ea79SLois Curfman McInnes 
5848965ea79SLois Curfman McInnes     if (!rank) {
585878740d9SKris Buschelman       ierr = MatCreate(mat->comm,M,N,M,N,&A);CHKERRQ(ierr);
5863a40ed3dSBarry Smith     } else {
587878740d9SKris Buschelman       ierr = MatCreate(mat->comm,0,0,M,N,&A);CHKERRQ(ierr);
5888965ea79SLois Curfman McInnes     }
589be5d1d56SKris Buschelman     /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */
590878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
591878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
592*52e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
5938965ea79SLois Curfman McInnes 
59439ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
59539ddd567SLois Curfman McInnes        but it's quick for now */
59651022da4SBarry Smith     A->insertmode = INSERT_VALUES;
597273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
5988965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
599ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
600ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
601ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
60239ddd567SLois Curfman McInnes       row++;
6038965ea79SLois Curfman McInnes     }
6048965ea79SLois Curfman McInnes 
6056d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6066d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
607b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
608b9b97703SBarry Smith     if (!rank) {
6096831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
6108965ea79SLois Curfman McInnes     }
611b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
612b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6138965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
6148965ea79SLois Curfman McInnes   }
6153a40ed3dSBarry Smith   PetscFunctionReturn(0);
6168965ea79SLois Curfman McInnes }
6178965ea79SLois Curfman McInnes 
6184a2ae208SSatish Balay #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
620dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
6218965ea79SLois Curfman McInnes {
622dfbe8321SBarry Smith   PetscErrorCode ierr;
62332077d6dSBarry Smith   PetscTruth iascii,isbinary,isdraw,issocket;
6248965ea79SLois Curfman McInnes 
625433994e6SBarry Smith   PetscFunctionBegin;
6260f5bd95cSBarry Smith 
62732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
628fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
629b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
630fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6310f5bd95cSBarry Smith 
63232077d6dSBarry Smith   if (iascii || issocket || isdraw) {
633f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6340f5bd95cSBarry Smith   } else if (isbinary) {
6353a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6365cd90555SBarry Smith   } else {
637958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6388965ea79SLois Curfman McInnes   }
6393a40ed3dSBarry Smith   PetscFunctionReturn(0);
6408965ea79SLois Curfman McInnes }
6418965ea79SLois Curfman McInnes 
6424a2ae208SSatish Balay #undef __FUNCT__
6434a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
644dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6458965ea79SLois Curfman McInnes {
6463501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6473501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
648dfbe8321SBarry Smith   PetscErrorCode ierr;
649329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6508965ea79SLois Curfman McInnes 
6513a40ed3dSBarry Smith   PetscFunctionBegin;
652273d9f13SBarry Smith   info->rows_global    = (double)A->M;
653273d9f13SBarry Smith   info->columns_global = (double)A->N;
654273d9f13SBarry Smith   info->rows_local     = (double)A->m;
655273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6564e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6574e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6584e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6594e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6608965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6614e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6624e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6634e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6644e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6654e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6668965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
667d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
6684e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6694e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6704e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6714e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6724e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6738965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
674d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
6754e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6764e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6774e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6784e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6794e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6808965ea79SLois Curfman McInnes   }
6814e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6824e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6834e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6843a40ed3dSBarry Smith   PetscFunctionReturn(0);
6858965ea79SLois Curfman McInnes }
6868965ea79SLois Curfman McInnes 
6874a2ae208SSatish Balay #undef __FUNCT__
6884a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
689dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op)
6908965ea79SLois Curfman McInnes {
69139ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
692dfbe8321SBarry Smith   PetscErrorCode ierr;
6938965ea79SLois Curfman McInnes 
6943a40ed3dSBarry Smith   PetscFunctionBegin;
69512c028f9SKris Buschelman   switch (op) {
69612c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
69712c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
69812c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
69912c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
70012c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
70112c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
702273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
70312c028f9SKris Buschelman     break;
70412c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
705273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
706273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
70712c028f9SKris Buschelman     break;
70812c028f9SKris Buschelman   case MAT_ROWS_SORTED:
70912c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
71012c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
71112c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
712b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
71312c028f9SKris Buschelman     break;
71412c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
715273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
716273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
71712c028f9SKris Buschelman     break;
71812c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
719273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
72012c028f9SKris Buschelman     break;
72112c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
72229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
72377e54ba9SKris Buschelman   case MAT_SYMMETRIC:
72477e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
7259a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
7269a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
7279a4540c5SBarry Smith   case MAT_HERMITIAN:
7289a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
7299a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
7309a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
73177e54ba9SKris Buschelman     break;
73212c028f9SKris Buschelman   default:
73329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
7343a40ed3dSBarry Smith   }
7353a40ed3dSBarry Smith   PetscFunctionReturn(0);
7368965ea79SLois Curfman McInnes }
7378965ea79SLois Curfman McInnes 
7388965ea79SLois Curfman McInnes 
7394a2ae208SSatish Balay #undef __FUNCT__
7404a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
741dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7425b2fa520SLois Curfman McInnes {
7435b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7445b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
74587828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
746dfbe8321SBarry Smith   PetscErrorCode ierr;
74713f74950SBarry Smith   PetscInt          i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7485b2fa520SLois Curfman McInnes 
7495b2fa520SLois Curfman McInnes   PetscFunctionBegin;
75072d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7515b2fa520SLois Curfman McInnes   if (ll) {
75272d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
75329bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
7541ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7555b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7565b2fa520SLois Curfman McInnes       x = l[i];
7575b2fa520SLois Curfman McInnes       v = mat->v + i;
7585b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7595b2fa520SLois Curfman McInnes     }
7601ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
761b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7625b2fa520SLois Curfman McInnes   }
7635b2fa520SLois Curfman McInnes   if (rr) {
76472d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
76529bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7665b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7675b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7681ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
7695b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7705b2fa520SLois Curfman McInnes       x = r[i];
7715b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7725b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7735b2fa520SLois Curfman McInnes     }
7741ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
775b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7765b2fa520SLois Curfman McInnes   }
7775b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7785b2fa520SLois Curfman McInnes }
7795b2fa520SLois Curfman McInnes 
7804a2ae208SSatish Balay #undef __FUNCT__
7814a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
782dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
783096963f5SLois Curfman McInnes {
7843501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7853501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
786dfbe8321SBarry Smith   PetscErrorCode ierr;
78713f74950SBarry Smith   PetscInt          i,j;
788329f5518SBarry Smith   PetscReal    sum = 0.0;
78987828ca2SBarry Smith   PetscScalar  *v = mat->v;
7903501a2bdSLois Curfman McInnes 
7913a40ed3dSBarry Smith   PetscFunctionBegin;
7923501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
793064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
7943501a2bdSLois Curfman McInnes   } else {
7953501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
796273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
797aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
798329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
7993501a2bdSLois Curfman McInnes #else
8003501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
8013501a2bdSLois Curfman McInnes #endif
8023501a2bdSLois Curfman McInnes       }
803064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
804064f8208SBarry Smith       *nrm = sqrt(*nrm);
805b0a32e0cSBarry Smith       PetscLogFlops(2*mdn->A->n*mdn->A->m);
8063a40ed3dSBarry Smith     } else if (type == NORM_1) {
807329f5518SBarry Smith       PetscReal *tmp,*tmp2;
808b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
809273d9f13SBarry Smith       tmp2 = tmp + A->N;
810273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
811064f8208SBarry Smith       *nrm = 0.0;
8123501a2bdSLois Curfman McInnes       v = mat->v;
813273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
814273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
81567e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8163501a2bdSLois Curfman McInnes         }
8173501a2bdSLois Curfman McInnes       }
818d7d1e502SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
819273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
820064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8213501a2bdSLois Curfman McInnes       }
822606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
823b0a32e0cSBarry Smith       PetscLogFlops(A->n*A->m);
8243a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
825329f5518SBarry Smith       PetscReal ntemp;
8263501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
827064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
8283a40ed3dSBarry Smith     } else {
82929bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8303501a2bdSLois Curfman McInnes     }
8313501a2bdSLois Curfman McInnes   }
8323a40ed3dSBarry Smith   PetscFunctionReturn(0);
8333501a2bdSLois Curfman McInnes }
8343501a2bdSLois Curfman McInnes 
8354a2ae208SSatish Balay #undef __FUNCT__
8364a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
837dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout)
8383501a2bdSLois Curfman McInnes {
8393501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8403501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8413501a2bdSLois Curfman McInnes   Mat          B;
84213f74950SBarry Smith   PetscInt          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8436849ba73SBarry Smith   PetscErrorCode ierr;
84413f74950SBarry Smith   PetscInt          j,i;
84587828ca2SBarry Smith   PetscScalar  *v;
8463501a2bdSLois Curfman McInnes 
8473a40ed3dSBarry Smith   PetscFunctionBegin;
8487c922b88SBarry Smith   if (!matout && M != N) {
84929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8507056b6fcSBarry Smith   }
851878740d9SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B);CHKERRQ(ierr);
852878740d9SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
853878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
8543501a2bdSLois Curfman McInnes 
855273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
85613f74950SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
8573501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8583501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8593501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8603501a2bdSLois Curfman McInnes     v   += m;
8613501a2bdSLois Curfman McInnes   }
862606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8636d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8646d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8657c922b88SBarry Smith   if (matout) {
8663501a2bdSLois Curfman McInnes     *matout = B;
8673501a2bdSLois Curfman McInnes   } else {
868273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8693501a2bdSLois Curfman McInnes   }
8703a40ed3dSBarry Smith   PetscFunctionReturn(0);
871096963f5SLois Curfman McInnes }
872096963f5SLois Curfman McInnes 
873d9eff348SSatish Balay #include "petscblaslapack.h"
8744a2ae208SSatish Balay #undef __FUNCT__
8754a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
876dfbe8321SBarry Smith PetscErrorCode MatScale_MPIDense(const PetscScalar *alpha,Mat inA)
87744cd7ae7SLois Curfman McInnes {
87844cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
87944cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
8804ce68768SBarry Smith   PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N;
88144cd7ae7SLois Curfman McInnes 
8823a40ed3dSBarry Smith   PetscFunctionBegin;
88371044d3cSBarry Smith   BLASscal_(&nz,(PetscScalar*)alpha,a->v,&one);
884b0a32e0cSBarry Smith   PetscLogFlops(nz);
8853a40ed3dSBarry Smith   PetscFunctionReturn(0);
88644cd7ae7SLois Curfman McInnes }
88744cd7ae7SLois Curfman McInnes 
8886849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
8898965ea79SLois Curfman McInnes 
8904a2ae208SSatish Balay #undef __FUNCT__
8914a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
892dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
893273d9f13SBarry Smith {
894dfbe8321SBarry Smith   PetscErrorCode ierr;
895273d9f13SBarry Smith 
896273d9f13SBarry Smith   PetscFunctionBegin;
897273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
898273d9f13SBarry Smith   PetscFunctionReturn(0);
899273d9f13SBarry Smith }
900273d9f13SBarry Smith 
9018965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
90209dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
90309dc0095SBarry Smith        MatGetRow_MPIDense,
90409dc0095SBarry Smith        MatRestoreRow_MPIDense,
90509dc0095SBarry Smith        MatMult_MPIDense,
90697304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
9077c922b88SBarry Smith        MatMultTranspose_MPIDense,
9087c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9098965ea79SLois Curfman McInnes        0,
91009dc0095SBarry Smith        0,
91109dc0095SBarry Smith        0,
91297304618SKris Buschelman /*10*/ 0,
91309dc0095SBarry Smith        0,
91409dc0095SBarry Smith        0,
91509dc0095SBarry Smith        0,
91609dc0095SBarry Smith        MatTranspose_MPIDense,
91797304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
91897304618SKris Buschelman        0,
91909dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9205b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
92109dc0095SBarry Smith        MatNorm_MPIDense,
92297304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
92309dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
92409dc0095SBarry Smith        0,
92509dc0095SBarry Smith        MatSetOption_MPIDense,
92609dc0095SBarry Smith        MatZeroEntries_MPIDense,
92797304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
92809dc0095SBarry Smith        0,
92909dc0095SBarry Smith        0,
93009dc0095SBarry Smith        0,
93109dc0095SBarry Smith        0,
93297304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
933273d9f13SBarry Smith        0,
93409dc0095SBarry Smith        0,
93509dc0095SBarry Smith        MatGetArray_MPIDense,
93609dc0095SBarry Smith        MatRestoreArray_MPIDense,
93797304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
93809dc0095SBarry Smith        0,
93909dc0095SBarry Smith        0,
94009dc0095SBarry Smith        0,
94109dc0095SBarry Smith        0,
94297304618SKris Buschelman /*40*/ 0,
9432ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
94409dc0095SBarry Smith        0,
94509dc0095SBarry Smith        MatGetValues_MPIDense,
94609dc0095SBarry Smith        0,
94797304618SKris Buschelman /*45*/ 0,
94809dc0095SBarry Smith        MatScale_MPIDense,
94909dc0095SBarry Smith        0,
95009dc0095SBarry Smith        0,
95109dc0095SBarry Smith        0,
952521d7252SBarry Smith /*50*/ 0,
95309dc0095SBarry Smith        0,
95409dc0095SBarry Smith        0,
95509dc0095SBarry Smith        0,
95609dc0095SBarry Smith        0,
95797304618SKris Buschelman /*55*/ 0,
95809dc0095SBarry Smith        0,
95909dc0095SBarry Smith        0,
96009dc0095SBarry Smith        0,
96109dc0095SBarry Smith        0,
96297304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
963b9b97703SBarry Smith        MatDestroy_MPIDense,
964b9b97703SBarry Smith        MatView_MPIDense,
96597304618SKris Buschelman        MatGetPetscMaps_Petsc,
96697304618SKris Buschelman        0,
96797304618SKris Buschelman /*65*/ 0,
96897304618SKris Buschelman        0,
96997304618SKris Buschelman        0,
97097304618SKris Buschelman        0,
97197304618SKris Buschelman        0,
97297304618SKris Buschelman /*70*/ 0,
97397304618SKris Buschelman        0,
97497304618SKris Buschelman        0,
97597304618SKris Buschelman        0,
97697304618SKris Buschelman        0,
97797304618SKris Buschelman /*75*/ 0,
97897304618SKris Buschelman        0,
97997304618SKris Buschelman        0,
98097304618SKris Buschelman        0,
98197304618SKris Buschelman        0,
98297304618SKris Buschelman /*80*/ 0,
98397304618SKris Buschelman        0,
98497304618SKris Buschelman        0,
98597304618SKris Buschelman        0,
986865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense,
987865e5f61SKris Buschelman        0,
988865e5f61SKris Buschelman        0,
989865e5f61SKris Buschelman        0,
990865e5f61SKris Buschelman        0,
991865e5f61SKris Buschelman        0,
992865e5f61SKris Buschelman /*90*/ 0,
993865e5f61SKris Buschelman        0,
994865e5f61SKris Buschelman        0,
995865e5f61SKris Buschelman        0,
996865e5f61SKris Buschelman        0,
997865e5f61SKris Buschelman /*95*/ 0,
998865e5f61SKris Buschelman        0,
999865e5f61SKris Buschelman        0,
1000865e5f61SKris Buschelman        0};
10018965ea79SLois Curfman McInnes 
1002273d9f13SBarry Smith EXTERN_C_BEGIN
10034a2ae208SSatish Balay #undef __FUNCT__
1004a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
1005dfbe8321SBarry Smith PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1006a23d5eceSKris Buschelman {
1007a23d5eceSKris Buschelman   Mat_MPIDense *a;
1008dfbe8321SBarry Smith   PetscErrorCode ierr;
1009a23d5eceSKris Buschelman 
1010a23d5eceSKris Buschelman   PetscFunctionBegin;
1011a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1012a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1013a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1014a23d5eceSKris Buschelman 
1015a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
10165c5985e7SKris Buschelman   ierr = MatCreate(PETSC_COMM_SELF,mat->m,mat->N,mat->m,mat->N,&a->A);CHKERRQ(ierr);
10175c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
10185c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
1019*52e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1020a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1021a23d5eceSKris Buschelman }
1022a23d5eceSKris Buschelman EXTERN_C_END
1023a23d5eceSKris Buschelman 
10240bad9183SKris Buschelman /*MC
1025fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
10260bad9183SKris Buschelman 
10270bad9183SKris Buschelman    Options Database Keys:
10280bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
10290bad9183SKris Buschelman 
10300bad9183SKris Buschelman   Level: beginner
10310bad9183SKris Buschelman 
10320bad9183SKris Buschelman .seealso: MatCreateMPIDense
10330bad9183SKris Buschelman M*/
10340bad9183SKris Buschelman 
1035a23d5eceSKris Buschelman EXTERN_C_BEGIN
1036a23d5eceSKris Buschelman #undef __FUNCT__
10374a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1038dfbe8321SBarry Smith PetscErrorCode MatCreate_MPIDense(Mat mat)
1039273d9f13SBarry Smith {
1040273d9f13SBarry Smith   Mat_MPIDense *a;
1041dfbe8321SBarry Smith   PetscErrorCode ierr;
104213f74950SBarry Smith   PetscInt          i;
1043273d9f13SBarry Smith 
1044273d9f13SBarry Smith   PetscFunctionBegin;
1045b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1046b0a32e0cSBarry Smith   mat->data         = (void*)a;
1047273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1048273d9f13SBarry Smith   mat->factor       = 0;
1049273d9f13SBarry Smith   mat->mapping      = 0;
1050273d9f13SBarry Smith 
1051273d9f13SBarry Smith   a->factor       = 0;
1052273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1053273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
1054273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
1055273d9f13SBarry Smith 
1056273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
1057273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
1058273d9f13SBarry Smith   a->nvec = mat->n;
1059273d9f13SBarry Smith 
1060273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1061273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
10628a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
10638a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
1064273d9f13SBarry Smith 
1065273d9f13SBarry Smith   /* build local table of row and column ownerships */
106613f74950SBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr);
1067273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
1068*52e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr);
106913f74950SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr);
1070273d9f13SBarry Smith   a->rowners[0] = 0;
1071273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1072273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
1073273d9f13SBarry Smith   }
1074273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
1075273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
107613f74950SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr);
1077273d9f13SBarry Smith   a->cowners[0] = 0;
1078273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1079273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1080273d9f13SBarry Smith   }
1081273d9f13SBarry Smith 
1082273d9f13SBarry Smith   /* build cache for off array entries formed */
1083273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1084273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1085273d9f13SBarry Smith 
1086273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1087273d9f13SBarry Smith   a->lvec        = 0;
1088273d9f13SBarry Smith   a->Mvctx       = 0;
1089273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1090273d9f13SBarry Smith 
1091273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1092273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1093273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1094a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1095a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1096a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1097273d9f13SBarry Smith   PetscFunctionReturn(0);
1098273d9f13SBarry Smith }
1099273d9f13SBarry Smith EXTERN_C_END
1100273d9f13SBarry Smith 
1101209238afSKris Buschelman /*MC
1102002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1103209238afSKris Buschelman 
1104209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1105209238afSKris Buschelman    and MATMPIDENSE otherwise.
1106209238afSKris Buschelman 
1107209238afSKris Buschelman    Options Database Keys:
1108209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1109209238afSKris Buschelman 
1110209238afSKris Buschelman   Level: beginner
1111209238afSKris Buschelman 
1112209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1113209238afSKris Buschelman M*/
1114209238afSKris Buschelman 
1115209238afSKris Buschelman EXTERN_C_BEGIN
1116209238afSKris Buschelman #undef __FUNCT__
1117209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1118dfbe8321SBarry Smith PetscErrorCode MatCreate_Dense(Mat A)
1119dfbe8321SBarry Smith {
11206849ba73SBarry Smith   PetscErrorCode ierr;
112113f74950SBarry Smith   PetscMPIInt    size;
1122209238afSKris Buschelman 
1123209238afSKris Buschelman   PetscFunctionBegin;
1124209238afSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr);
1125209238afSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1126209238afSKris Buschelman   if (size == 1) {
1127209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1128209238afSKris Buschelman   } else {
1129209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1130209238afSKris Buschelman   }
1131209238afSKris Buschelman   PetscFunctionReturn(0);
1132209238afSKris Buschelman }
1133209238afSKris Buschelman EXTERN_C_END
1134209238afSKris Buschelman 
11354a2ae208SSatish Balay #undef __FUNCT__
11364a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1137273d9f13SBarry Smith /*@C
1138273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1139273d9f13SBarry Smith 
1140273d9f13SBarry Smith    Not collective
1141273d9f13SBarry Smith 
1142273d9f13SBarry Smith    Input Parameters:
1143273d9f13SBarry Smith .  A - the matrix
1144273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1145273d9f13SBarry Smith    to control all matrix memory allocation.
1146273d9f13SBarry Smith 
1147273d9f13SBarry Smith    Notes:
1148273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1149273d9f13SBarry Smith    storage by columns.
1150273d9f13SBarry Smith 
1151273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1152273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1153273d9f13SBarry Smith    set data=PETSC_NULL.
1154273d9f13SBarry Smith 
1155273d9f13SBarry Smith    Level: intermediate
1156273d9f13SBarry Smith 
1157273d9f13SBarry Smith .keywords: matrix,dense, parallel
1158273d9f13SBarry Smith 
1159273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1160273d9f13SBarry Smith @*/
1161dfbe8321SBarry Smith PetscErrorCode MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1162273d9f13SBarry Smith {
1163dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar *);
1164273d9f13SBarry Smith 
1165273d9f13SBarry Smith   PetscFunctionBegin;
1166565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1167a23d5eceSKris Buschelman   if (f) {
1168a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1169a23d5eceSKris Buschelman   }
1170273d9f13SBarry Smith   PetscFunctionReturn(0);
1171273d9f13SBarry Smith }
1172273d9f13SBarry Smith 
11734a2ae208SSatish Balay #undef __FUNCT__
11744a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
11758965ea79SLois Curfman McInnes /*@C
117639ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
11778965ea79SLois Curfman McInnes 
1178db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1179db81eaa0SLois Curfman McInnes 
11808965ea79SLois Curfman McInnes    Input Parameters:
1181db81eaa0SLois Curfman McInnes +  comm - MPI communicator
11828965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1183db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
11848965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1185db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
11867f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1187dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
11888965ea79SLois Curfman McInnes 
11898965ea79SLois Curfman McInnes    Output Parameter:
1190477f1c0bSLois Curfman McInnes .  A - the matrix
11918965ea79SLois Curfman McInnes 
1192b259b22eSLois Curfman McInnes    Notes:
119339ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
119439ddd567SLois Curfman McInnes    storage by columns.
11958965ea79SLois Curfman McInnes 
119618f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
119718f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
11987f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
119918f449edSLois Curfman McInnes 
12008965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
12018965ea79SLois Curfman McInnes    (possibly both).
12028965ea79SLois Curfman McInnes 
1203027ccd11SLois Curfman McInnes    Level: intermediate
1204027ccd11SLois Curfman McInnes 
120539ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
12068965ea79SLois Curfman McInnes 
120739ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
12088965ea79SLois Curfman McInnes @*/
120913f74950SBarry Smith PetscErrorCode MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
12108965ea79SLois Curfman McInnes {
12116849ba73SBarry Smith   PetscErrorCode ierr;
121213f74950SBarry Smith   PetscMPIInt    size;
12138965ea79SLois Curfman McInnes 
12143a40ed3dSBarry Smith   PetscFunctionBegin;
1215273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1216273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1217273d9f13SBarry Smith   if (size > 1) {
1218273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1219273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1220273d9f13SBarry Smith   } else {
1221273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1222273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
12238c469469SLois Curfman McInnes   }
12243a40ed3dSBarry Smith   PetscFunctionReturn(0);
12258965ea79SLois Curfman McInnes }
12268965ea79SLois Curfman McInnes 
12274a2ae208SSatish Balay #undef __FUNCT__
12284a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
12296849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12308965ea79SLois Curfman McInnes {
12318965ea79SLois Curfman McInnes   Mat          mat;
12323501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
1233dfbe8321SBarry Smith   PetscErrorCode ierr;
12348965ea79SLois Curfman McInnes 
12353a40ed3dSBarry Smith   PetscFunctionBegin;
12368965ea79SLois Curfman McInnes   *newmat       = 0;
1237273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1238be5d1d56SKris Buschelman   ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr);
1239b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1240b0a32e0cSBarry Smith   mat->data         = (void*)a;
1241549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
12423501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1243c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1244273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
12458965ea79SLois Curfman McInnes 
12468965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
12478965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
12488965ea79SLois Curfman McInnes   a->size         = oldmat->size;
12498965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1250e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1251b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
12523782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
125313f74950SBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr);
1254*52e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr);
125513f74950SBarry Smith   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
12568798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
12578965ea79SLois Curfman McInnes 
1258329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
12595609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
1260*52e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
12618965ea79SLois Curfman McInnes   *newmat = mat;
12623a40ed3dSBarry Smith   PetscFunctionReturn(0);
12638965ea79SLois Curfman McInnes }
12648965ea79SLois Curfman McInnes 
1265e090d566SSatish Balay #include "petscsys.h"
12668965ea79SLois Curfman McInnes 
12674a2ae208SSatish Balay #undef __FUNCT__
12684a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
126913f74950SBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,const MatType type,Mat *newmat)
127090ace30eSBarry Smith {
12716849ba73SBarry Smith   PetscErrorCode ierr;
127232dcc486SBarry Smith   PetscMPIInt    rank,size;
127313f74950SBarry Smith   PetscInt            *rowners,i,m,nz,j;
127487828ca2SBarry Smith   PetscScalar    *array,*vals,*vals_ptr;
127590ace30eSBarry Smith   MPI_Status     status;
127690ace30eSBarry Smith 
12773a40ed3dSBarry Smith   PetscFunctionBegin;
1278d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1279d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
128090ace30eSBarry Smith 
128190ace30eSBarry Smith   /* determine ownership of all rows */
128290ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
128313f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
128413f74950SBarry Smith   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
128590ace30eSBarry Smith   rowners[0] = 0;
128690ace30eSBarry Smith   for (i=2; i<=size; i++) {
128790ace30eSBarry Smith     rowners[i] += rowners[i-1];
128890ace30eSBarry Smith   }
128990ace30eSBarry Smith 
1290878740d9SKris Buschelman   ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr);
1291be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1292878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
129390ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
129490ace30eSBarry Smith 
129590ace30eSBarry Smith   if (!rank) {
129687828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
129790ace30eSBarry Smith 
129890ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
12990752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
130090ace30eSBarry Smith 
130190ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
130290ace30eSBarry Smith     vals_ptr = vals;
130390ace30eSBarry Smith     for (i=0; i<m; i++) {
130490ace30eSBarry Smith       for (j=0; j<N; j++) {
130590ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
130690ace30eSBarry Smith       }
130790ace30eSBarry Smith     }
130890ace30eSBarry Smith 
130990ace30eSBarry Smith     /* read in other processors and ship out */
131090ace30eSBarry Smith     for (i=1; i<size; i++) {
131190ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
13120752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1313ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
131490ace30eSBarry Smith     }
13153a40ed3dSBarry Smith   } else {
131690ace30eSBarry Smith     /* receive numeric values */
131787828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
131890ace30eSBarry Smith 
131990ace30eSBarry Smith     /* receive message of values*/
1320ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
132190ace30eSBarry Smith 
132290ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
132390ace30eSBarry Smith     vals_ptr = vals;
132490ace30eSBarry Smith     for (i=0; i<m; i++) {
132590ace30eSBarry Smith       for (j=0; j<N; j++) {
132690ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
132790ace30eSBarry Smith       }
132890ace30eSBarry Smith     }
132990ace30eSBarry Smith   }
1330606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1331606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
13326d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13336d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13343a40ed3dSBarry Smith   PetscFunctionReturn(0);
133590ace30eSBarry Smith }
133690ace30eSBarry Smith 
13374a2ae208SSatish Balay #undef __FUNCT__
13384a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
1339dfbe8321SBarry Smith PetscErrorCode MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat)
13408965ea79SLois Curfman McInnes {
13418965ea79SLois Curfman McInnes   Mat            A;
134287828ca2SBarry Smith   PetscScalar    *vals,*svals;
134319bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
13448965ea79SLois Curfman McInnes   MPI_Status     status;
134513f74950SBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
134613f74950SBarry Smith   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
134713f74950SBarry Smith   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
134813f74950SBarry Smith   PetscInt       i,nz,j,rstart,rend;
134913f74950SBarry Smith   int            fd;
13506849ba73SBarry Smith   PetscErrorCode ierr;
13518965ea79SLois Curfman McInnes 
13523a40ed3dSBarry Smith   PetscFunctionBegin;
1353d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1354d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13558965ea79SLois Curfman McInnes   if (!rank) {
1356b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
13570752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1358552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
13598965ea79SLois Curfman McInnes   }
13608965ea79SLois Curfman McInnes 
136113f74950SBarry Smith   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
136290ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
136390ace30eSBarry Smith 
136490ace30eSBarry Smith   /*
136590ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
136690ace30eSBarry Smith   */
136790ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1368be5d1d56SKris Buschelman     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr);
13693a40ed3dSBarry Smith     PetscFunctionReturn(0);
137090ace30eSBarry Smith   }
137190ace30eSBarry Smith 
13728965ea79SLois Curfman McInnes   /* determine ownership of all rows */
13738965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
137413f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
1375ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
13768965ea79SLois Curfman McInnes   rowners[0] = 0;
13778965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
13788965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
13798965ea79SLois Curfman McInnes   }
13808965ea79SLois Curfman McInnes   rstart = rowners[rank];
13818965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
13828965ea79SLois Curfman McInnes 
13838965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
138413f74950SBarry Smith   ierr    = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr);
13858965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
13868965ea79SLois Curfman McInnes   if (!rank) {
138713f74950SBarry Smith     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
13880752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
138913f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
13908965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
13911466f1e1SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
1392606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1393ca161407SBarry Smith   } else {
13941466f1e1SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
13958965ea79SLois Curfman McInnes   }
13968965ea79SLois Curfman McInnes 
13978965ea79SLois Curfman McInnes   if (!rank) {
13988965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
139913f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
140013f74950SBarry Smith     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
14018965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14028965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
14038965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
14048965ea79SLois Curfman McInnes       }
14058965ea79SLois Curfman McInnes     }
1406606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
14078965ea79SLois Curfman McInnes 
14088965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
14098965ea79SLois Curfman McInnes     maxnz = 0;
14108965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14110452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
14128965ea79SLois Curfman McInnes     }
141313f74950SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
14148965ea79SLois Curfman McInnes 
14158965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
14168965ea79SLois Curfman McInnes     nz = procsnz[0];
141713f74950SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
14180752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
14198965ea79SLois Curfman McInnes 
14208965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
14218965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14228965ea79SLois Curfman McInnes       nz   = procsnz[i];
14230752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
142413f74950SBarry Smith       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
14258965ea79SLois Curfman McInnes     }
1426606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
14273a40ed3dSBarry Smith   } else {
14288965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
14298965ea79SLois Curfman McInnes     nz = 0;
14308965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14318965ea79SLois Curfman McInnes       nz += ourlens[i];
14328965ea79SLois Curfman McInnes     }
143313f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
14348965ea79SLois Curfman McInnes 
14358965ea79SLois Curfman McInnes     /* receive message of column indices*/
143613f74950SBarry Smith     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
143713f74950SBarry Smith     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
143829bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14398965ea79SLois Curfman McInnes   }
14408965ea79SLois Curfman McInnes 
14418965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
144213f74950SBarry Smith   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
14438965ea79SLois Curfman McInnes   jj = 0;
14448965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14458965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
14468965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
14478965ea79SLois Curfman McInnes       jj++;
14488965ea79SLois Curfman McInnes     }
14498965ea79SLois Curfman McInnes   }
14508965ea79SLois Curfman McInnes 
14518965ea79SLois Curfman McInnes   /* create our matrix */
14528965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14538965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
14548965ea79SLois Curfman McInnes   }
1455878740d9SKris Buschelman   ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr);
1456be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1457878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
14588965ea79SLois Curfman McInnes   A = *newmat;
14598965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14608965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
14618965ea79SLois Curfman McInnes   }
14628965ea79SLois Curfman McInnes 
14638965ea79SLois Curfman McInnes   if (!rank) {
146487828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14658965ea79SLois Curfman McInnes 
14668965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
14678965ea79SLois Curfman McInnes     nz = procsnz[0];
14680752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
14698965ea79SLois Curfman McInnes 
14708965ea79SLois Curfman McInnes     /* insert into matrix */
14718965ea79SLois Curfman McInnes     jj      = rstart;
14728965ea79SLois Curfman McInnes     smycols = mycols;
14738965ea79SLois Curfman McInnes     svals   = vals;
14748965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14758965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14768965ea79SLois Curfman McInnes       smycols += ourlens[i];
14778965ea79SLois Curfman McInnes       svals   += ourlens[i];
14788965ea79SLois Curfman McInnes       jj++;
14798965ea79SLois Curfman McInnes     }
14808965ea79SLois Curfman McInnes 
14818965ea79SLois Curfman McInnes     /* read in other processors and ship out */
14828965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14838965ea79SLois Curfman McInnes       nz   = procsnz[i];
14840752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1485ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
14868965ea79SLois Curfman McInnes     }
1487606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
14883a40ed3dSBarry Smith   } else {
14898965ea79SLois Curfman McInnes     /* receive numeric values */
149084d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14918965ea79SLois Curfman McInnes 
14928965ea79SLois Curfman McInnes     /* receive message of values*/
1493ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1494ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
149529bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14968965ea79SLois Curfman McInnes 
14978965ea79SLois Curfman McInnes     /* insert into matrix */
14988965ea79SLois Curfman McInnes     jj      = rstart;
14998965ea79SLois Curfman McInnes     smycols = mycols;
15008965ea79SLois Curfman McInnes     svals   = vals;
15018965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15028965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15038965ea79SLois Curfman McInnes       smycols += ourlens[i];
15048965ea79SLois Curfman McInnes       svals   += ourlens[i];
15058965ea79SLois Curfman McInnes       jj++;
15068965ea79SLois Curfman McInnes     }
15078965ea79SLois Curfman McInnes   }
1508606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1509606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1510606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1511606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
15128965ea79SLois Curfman McInnes 
15136d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15146d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15153a40ed3dSBarry Smith   PetscFunctionReturn(0);
15168965ea79SLois Curfman McInnes }
151790ace30eSBarry Smith 
151890ace30eSBarry Smith 
151990ace30eSBarry Smith 
152090ace30eSBarry Smith 
152190ace30eSBarry Smith 
1522