xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 89ac3759e191d93affd9fbb533c0e69bdb54a93d)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3ed3cc1f0SBarry Smith /*
4ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
5ed3cc1f0SBarry Smith */
6ed3cc1f0SBarry Smith 
7*89ac3759SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"    /*I   "petscmat.h"  I*/
88965ea79SLois Curfman McInnes 
9ba8c8a56SBarry Smith #undef __FUNCT__
10ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
11ab92ecdeSBarry Smith /*@
12ab92ecdeSBarry Smith 
13ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
14ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
15ab92ecdeSBarry Smith 
16ab92ecdeSBarry Smith     Input Parameter:
17ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
18ab92ecdeSBarry Smith 
19ab92ecdeSBarry Smith     Output Parameter:
20ab92ecdeSBarry Smith .      B - the inner matrix
21ab92ecdeSBarry Smith 
22ab92ecdeSBarry Smith @*/
23ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
24ab92ecdeSBarry Smith {
25ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
26ab92ecdeSBarry Smith   PetscErrorCode ierr;
27ab92ecdeSBarry Smith   PetscTruth     flg;
28ab92ecdeSBarry Smith   PetscMPIInt    size;
29ab92ecdeSBarry Smith 
30ab92ecdeSBarry Smith   PetscFunctionBegin;
31ab92ecdeSBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
32ab92ecdeSBarry Smith   if (!flg) {  /* this check sucks! */
33ab92ecdeSBarry Smith     ierr = PetscTypeCompare((PetscObject)A,MATDENSE,&flg);CHKERRQ(ierr);
34ab92ecdeSBarry Smith     if (flg) {
35ab92ecdeSBarry Smith       ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
36ab92ecdeSBarry Smith       if (size == 1) flg = PETSC_FALSE;
37ab92ecdeSBarry Smith     }
38ab92ecdeSBarry Smith   }
39ab92ecdeSBarry Smith   if (flg) {
40ab92ecdeSBarry Smith     *B = mat->A;
41ab92ecdeSBarry Smith   } else {
42ab92ecdeSBarry Smith     *B = A;
43ab92ecdeSBarry Smith   }
44ab92ecdeSBarry Smith   PetscFunctionReturn(0);
45ab92ecdeSBarry Smith }
46ab92ecdeSBarry Smith 
47ab92ecdeSBarry Smith #undef __FUNCT__
48ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
49ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
50ba8c8a56SBarry Smith {
51ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
52ba8c8a56SBarry Smith   PetscErrorCode ierr;
53ba8c8a56SBarry Smith   PetscInt       lrow,rstart = mat->rstart,rend = mat->rend;
54ba8c8a56SBarry Smith 
55ba8c8a56SBarry Smith   PetscFunctionBegin;
56ba8c8a56SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
57ba8c8a56SBarry Smith   lrow = row - rstart;
58ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
59ba8c8a56SBarry Smith   PetscFunctionReturn(0);
60ba8c8a56SBarry Smith }
61ba8c8a56SBarry Smith 
62ba8c8a56SBarry Smith #undef __FUNCT__
63ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
64ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
65ba8c8a56SBarry Smith {
66ba8c8a56SBarry Smith   PetscErrorCode ierr;
67ba8c8a56SBarry Smith 
68ba8c8a56SBarry Smith   PetscFunctionBegin;
69ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
70ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
71ba8c8a56SBarry Smith   PetscFunctionReturn(0);
72ba8c8a56SBarry Smith }
73ba8c8a56SBarry Smith 
740de54da6SSatish Balay EXTERN_C_BEGIN
754a2ae208SSatish Balay #undef __FUNCT__
764a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
77be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
780de54da6SSatish Balay {
790de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
806849ba73SBarry Smith   PetscErrorCode ierr;
8113f74950SBarry Smith   PetscInt       m = A->m,rstart = mdn->rstart;
8287828ca2SBarry Smith   PetscScalar    *array;
830de54da6SSatish Balay   MPI_Comm       comm;
840de54da6SSatish Balay 
850de54da6SSatish Balay   PetscFunctionBegin;
86273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
870de54da6SSatish Balay 
880de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
890de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
900de54da6SSatish Balay 
910de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
920de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
93f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,B);CHKERRQ(ierr);
94f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr);
955c5985e7SKris Buschelman   ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr);
965c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
970de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
980de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
990de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1000de54da6SSatish Balay 
1010de54da6SSatish Balay   *iscopy = PETSC_TRUE;
1020de54da6SSatish Balay   PetscFunctionReturn(0);
1030de54da6SSatish Balay }
1040de54da6SSatish Balay EXTERN_C_END
1050de54da6SSatish Balay 
1064a2ae208SSatish Balay #undef __FUNCT__
1074a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
10813f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1098965ea79SLois Curfman McInnes {
11039b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
111dfbe8321SBarry Smith   PetscErrorCode ierr;
11213f74950SBarry Smith   PetscInt       i,j,rstart = A->rstart,rend = A->rend,row;
113273d9f13SBarry Smith   PetscTruth     roworiented = A->roworiented;
1148965ea79SLois Curfman McInnes 
1153a40ed3dSBarry Smith   PetscFunctionBegin;
1168965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1175ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
118273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1198965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1208965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
12139b7565bSBarry Smith       if (roworiented) {
12239b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1233a40ed3dSBarry Smith       } else {
1248965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1255ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
126273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
12739b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
12839b7565bSBarry Smith         }
1298965ea79SLois Curfman McInnes       }
1303a40ed3dSBarry Smith     } else {
1313782ba37SSatish Balay       if (!A->donotstash) {
13239b7565bSBarry Smith         if (roworiented) {
1338798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
134d36fbae8SSatish Balay         } else {
1358798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
13639b7565bSBarry Smith         }
137b49de8d1SLois Curfman McInnes       }
138b49de8d1SLois Curfman McInnes     }
1393782ba37SSatish Balay   }
1403a40ed3dSBarry Smith   PetscFunctionReturn(0);
141b49de8d1SLois Curfman McInnes }
142b49de8d1SLois Curfman McInnes 
1434a2ae208SSatish Balay #undef __FUNCT__
1444a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
14513f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
146b49de8d1SLois Curfman McInnes {
147b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
148dfbe8321SBarry Smith   PetscErrorCode ierr;
14913f74950SBarry Smith   PetscInt       i,j,rstart = mdn->rstart,rend = mdn->rend,row;
150b49de8d1SLois Curfman McInnes 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
152b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
15329bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
154273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
155b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
156b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
157b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
15829bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
159273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
16029bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
161a8c6a408SBarry Smith         }
162b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
163b49de8d1SLois Curfman McInnes       }
164a8c6a408SBarry Smith     } else {
16529bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
1668965ea79SLois Curfman McInnes     }
1678965ea79SLois Curfman McInnes   }
1683a40ed3dSBarry Smith   PetscFunctionReturn(0);
1698965ea79SLois Curfman McInnes }
1708965ea79SLois Curfman McInnes 
1714a2ae208SSatish Balay #undef __FUNCT__
1724a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
173dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
174ff14e315SSatish Balay {
175ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
176dfbe8321SBarry Smith   PetscErrorCode ierr;
177ff14e315SSatish Balay 
1783a40ed3dSBarry Smith   PetscFunctionBegin;
179ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1803a40ed3dSBarry Smith   PetscFunctionReturn(0);
181ff14e315SSatish Balay }
182ff14e315SSatish Balay 
1834a2ae208SSatish Balay #undef __FUNCT__
1844a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
18513f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
186ca3fa75bSLois Curfman McInnes {
187ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data,*newmatd;
188ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
1896849ba73SBarry Smith   PetscErrorCode ierr;
19013f74950SBarry Smith   PetscInt       i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
19187828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
192ca3fa75bSLois Curfman McInnes   Mat            newmat;
193ca3fa75bSLois Curfman McInnes 
194ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
195ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
196ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
197b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
198b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
199ca3fa75bSLois Curfman McInnes 
200ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
2017eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
2027eba5e9cSLois Curfman McInnes      original matrix! */
203ca3fa75bSLois Curfman McInnes 
204ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
2057eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
206ca3fa75bSLois Curfman McInnes 
207ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
208ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
20929bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2107eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
211ca3fa75bSLois Curfman McInnes     newmat = *B;
212ca3fa75bSLois Curfman McInnes   } else {
213ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
214f69a0ea3SMatthew Knepley     ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr);
215f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr);
216878740d9SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
217878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
218ca3fa75bSLois Curfman McInnes   }
219ca3fa75bSLois Curfman McInnes 
220ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
221ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
222ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
223ca3fa75bSLois Curfman McInnes 
224ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
225ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
226ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2277eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
228ca3fa75bSLois Curfman McInnes     }
229ca3fa75bSLois Curfman McInnes   }
230ca3fa75bSLois Curfman McInnes 
231ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
232ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
233ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
234ca3fa75bSLois Curfman McInnes 
235ca3fa75bSLois Curfman McInnes   /* Free work space */
236ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
237ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
238ca3fa75bSLois Curfman McInnes   *B = newmat;
239ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
240ca3fa75bSLois Curfman McInnes }
241ca3fa75bSLois Curfman McInnes 
2424a2ae208SSatish Balay #undef __FUNCT__
2434a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
244dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
245ff14e315SSatish Balay {
2463a40ed3dSBarry Smith   PetscFunctionBegin;
2473a40ed3dSBarry Smith   PetscFunctionReturn(0);
248ff14e315SSatish Balay }
249ff14e315SSatish Balay 
2504a2ae208SSatish Balay #undef __FUNCT__
2514a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
252dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2538965ea79SLois Curfman McInnes {
25439ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
2558965ea79SLois Curfman McInnes   MPI_Comm       comm = mat->comm;
256dfbe8321SBarry Smith   PetscErrorCode ierr;
25713f74950SBarry Smith   PetscInt       nstash,reallocs;
2588965ea79SLois Curfman McInnes   InsertMode     addv;
2598965ea79SLois Curfman McInnes 
2603a40ed3dSBarry Smith   PetscFunctionBegin;
2618965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
262ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
2637056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
26429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
2658965ea79SLois Curfman McInnes   }
266e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2678965ea79SLois Curfman McInnes 
2688798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
2698798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
27063ba0a88SBarry Smith   ierr = PetscLogInfo((mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs));CHKERRQ(ierr);
2713a40ed3dSBarry Smith   PetscFunctionReturn(0);
2728965ea79SLois Curfman McInnes }
2738965ea79SLois Curfman McInnes 
2744a2ae208SSatish Balay #undef __FUNCT__
2754a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
276dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2778965ea79SLois Curfman McInnes {
27839ddd567SLois Curfman McInnes   Mat_MPIDense    *mdn=(Mat_MPIDense*)mat->data;
2796849ba73SBarry Smith   PetscErrorCode  ierr;
28013f74950SBarry Smith   PetscInt        i,*row,*col,flg,j,rstart,ncols;
28113f74950SBarry Smith   PetscMPIInt     n;
28287828ca2SBarry Smith   PetscScalar     *val;
283e0fa3b82SLois Curfman McInnes   InsertMode      addv=mat->insertmode;
2848965ea79SLois Curfman McInnes 
2853a40ed3dSBarry Smith   PetscFunctionBegin;
2868965ea79SLois Curfman McInnes   /*  wait on receives */
2877ef1d9bdSSatish Balay   while (1) {
2888798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2897ef1d9bdSSatish Balay     if (!flg) break;
2908965ea79SLois Curfman McInnes 
2917ef1d9bdSSatish Balay     for (i=0; i<n;) {
2927ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2937ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2947ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2957ef1d9bdSSatish Balay       else       ncols = n-i;
2967ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2977ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2987ef1d9bdSSatish Balay       i = j;
2998965ea79SLois Curfman McInnes     }
3007ef1d9bdSSatish Balay   }
3018798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
3028965ea79SLois Curfman McInnes 
30339ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
30439ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3058965ea79SLois Curfman McInnes 
3066d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
30739ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3088965ea79SLois Curfman McInnes   }
3093a40ed3dSBarry Smith   PetscFunctionReturn(0);
3108965ea79SLois Curfman McInnes }
3118965ea79SLois Curfman McInnes 
3124a2ae208SSatish Balay #undef __FUNCT__
3134a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
314dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3158965ea79SLois Curfman McInnes {
316dfbe8321SBarry Smith   PetscErrorCode ierr;
31739ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3183a40ed3dSBarry Smith 
3193a40ed3dSBarry Smith   PetscFunctionBegin;
3203a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3213a40ed3dSBarry Smith   PetscFunctionReturn(0);
3228965ea79SLois Curfman McInnes }
3238965ea79SLois Curfman McInnes 
3248965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3258965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3268965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3273501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3288965ea79SLois Curfman McInnes    routine.
3298965ea79SLois Curfman McInnes */
3304a2ae208SSatish Balay #undef __FUNCT__
3314a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
332f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
3338965ea79SLois Curfman McInnes {
33439ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3356849ba73SBarry Smith   PetscErrorCode ierr;
336f4df32b1SMatthew Knepley   PetscInt       i,*owners = l->rowners;
33713f74950SBarry Smith   PetscInt       *nprocs,j,idx,nsends;
33813f74950SBarry Smith   PetscInt       nmax,*svalues,*starts,*owner,nrecvs;
33913f74950SBarry Smith   PetscInt       *rvalues,tag = A->tag,count,base,slen,*source;
34013f74950SBarry Smith   PetscInt       *lens,*lrows,*values;
34113f74950SBarry Smith   PetscMPIInt    n,imdex,rank = l->rank,size = l->size;
3428965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
3438965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
3448965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
34535d8aa7fSBarry Smith   PetscTruth     found;
3468965ea79SLois Curfman McInnes 
3473a40ed3dSBarry Smith   PetscFunctionBegin;
3488965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
34913f74950SBarry Smith   ierr  = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
35013f74950SBarry Smith   ierr  = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
35113f74950SBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/
3528965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3538965ea79SLois Curfman McInnes     idx = rows[i];
35435d8aa7fSBarry Smith     found = PETSC_FALSE;
3558965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3568965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
357c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3588965ea79SLois Curfman McInnes       }
3598965ea79SLois Curfman McInnes     }
36029bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3618965ea79SLois Curfman McInnes   }
362c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3638965ea79SLois Curfman McInnes 
3648965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
365c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3668965ea79SLois Curfman McInnes 
3678965ea79SLois Curfman McInnes   /* post receives:   */
36813f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
369b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3708965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
37113f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3728965ea79SLois Curfman McInnes   }
3738965ea79SLois Curfman McInnes 
3748965ea79SLois Curfman McInnes   /* do sends:
3758965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3768965ea79SLois Curfman McInnes          the ith processor
3778965ea79SLois Curfman McInnes   */
37813f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
379b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
38013f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
3818965ea79SLois Curfman McInnes   starts[0]  = 0;
382c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3838965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3848965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3858965ea79SLois Curfman McInnes   }
3868965ea79SLois Curfman McInnes 
3878965ea79SLois Curfman McInnes   starts[0] = 0;
388c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3898965ea79SLois Curfman McInnes   count = 0;
3908965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
391c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
39213f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3938965ea79SLois Curfman McInnes     }
3948965ea79SLois Curfman McInnes   }
395606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3968965ea79SLois Curfman McInnes 
3978965ea79SLois Curfman McInnes   base = owners[rank];
3988965ea79SLois Curfman McInnes 
3998965ea79SLois Curfman McInnes   /*  wait on receives */
40013f74950SBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
4018965ea79SLois Curfman McInnes   source = lens + nrecvs;
4028965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
4038965ea79SLois Curfman McInnes   while (count) {
404ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4058965ea79SLois Curfman McInnes     /* unpack receives into our local space */
40613f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4078965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
4088965ea79SLois Curfman McInnes     lens[imdex]    = n;
4098965ea79SLois Curfman McInnes     slen += n;
4108965ea79SLois Curfman McInnes     count--;
4118965ea79SLois Curfman McInnes   }
412606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4138965ea79SLois Curfman McInnes 
4148965ea79SLois Curfman McInnes   /* move the data into the send scatter */
41513f74950SBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
4168965ea79SLois Curfman McInnes   count = 0;
4178965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4188965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4198965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4208965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4218965ea79SLois Curfman McInnes     }
4228965ea79SLois Curfman McInnes   }
423606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
424606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
425606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
426606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
4278965ea79SLois Curfman McInnes 
4288965ea79SLois Curfman McInnes   /* actually zap the local rows */
429f4df32b1SMatthew Knepley   ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr);
430606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4318965ea79SLois Curfman McInnes 
4328965ea79SLois Curfman McInnes   /* wait on sends */
4338965ea79SLois Curfman McInnes   if (nsends) {
434b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
435ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
436606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4378965ea79SLois Curfman McInnes   }
438606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
439606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4408965ea79SLois Curfman McInnes 
4413a40ed3dSBarry Smith   PetscFunctionReturn(0);
4428965ea79SLois Curfman McInnes }
4438965ea79SLois Curfman McInnes 
4444a2ae208SSatish Balay #undef __FUNCT__
4454a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
446dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4478965ea79SLois Curfman McInnes {
44839ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
449dfbe8321SBarry Smith   PetscErrorCode ierr;
450c456f294SBarry Smith 
4513a40ed3dSBarry Smith   PetscFunctionBegin;
45243a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
45343a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
45444cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4553a40ed3dSBarry Smith   PetscFunctionReturn(0);
4568965ea79SLois Curfman McInnes }
4578965ea79SLois Curfman McInnes 
4584a2ae208SSatish Balay #undef __FUNCT__
4594a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
460dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4618965ea79SLois Curfman McInnes {
46239ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
463dfbe8321SBarry Smith   PetscErrorCode ierr;
464c456f294SBarry Smith 
4653a40ed3dSBarry Smith   PetscFunctionBegin;
46643a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
46743a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
46844cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4693a40ed3dSBarry Smith   PetscFunctionReturn(0);
4708965ea79SLois Curfman McInnes }
4718965ea79SLois Curfman McInnes 
4724a2ae208SSatish Balay #undef __FUNCT__
4734a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
474dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
475096963f5SLois Curfman McInnes {
476096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
477dfbe8321SBarry Smith   PetscErrorCode ierr;
47887828ca2SBarry Smith   PetscScalar    zero = 0.0;
479096963f5SLois Curfman McInnes 
4803a40ed3dSBarry Smith   PetscFunctionBegin;
4812dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
4827c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
483537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
484537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4853a40ed3dSBarry Smith   PetscFunctionReturn(0);
486096963f5SLois Curfman McInnes }
487096963f5SLois Curfman McInnes 
4884a2ae208SSatish Balay #undef __FUNCT__
4894a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
490dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
491096963f5SLois Curfman McInnes {
492096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
493dfbe8321SBarry Smith   PetscErrorCode ierr;
494096963f5SLois Curfman McInnes 
4953a40ed3dSBarry Smith   PetscFunctionBegin;
4963501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4977c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
498537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
499537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
5003a40ed3dSBarry Smith   PetscFunctionReturn(0);
501096963f5SLois Curfman McInnes }
502096963f5SLois Curfman McInnes 
5034a2ae208SSatish Balay #undef __FUNCT__
5044a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
505dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5068965ea79SLois Curfman McInnes {
50739ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
508096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
509dfbe8321SBarry Smith   PetscErrorCode ierr;
51013f74950SBarry Smith   PetscInt       len,i,n,m = A->m,radd;
51187828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
512ed3cc1f0SBarry Smith 
5133a40ed3dSBarry Smith   PetscFunctionBegin;
5142dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5151ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
516096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
517273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
518273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
5197ddc982cSLois Curfman McInnes   radd = a->rstart*m;
52044cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
521096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
522096963f5SLois Curfman McInnes   }
5231ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5243a40ed3dSBarry Smith   PetscFunctionReturn(0);
5258965ea79SLois Curfman McInnes }
5268965ea79SLois Curfman McInnes 
5274a2ae208SSatish Balay #undef __FUNCT__
5284a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
529dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5308965ea79SLois Curfman McInnes {
5313501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
532dfbe8321SBarry Smith   PetscErrorCode ierr;
533ed3cc1f0SBarry Smith 
5343a40ed3dSBarry Smith   PetscFunctionBegin;
53594d884c6SBarry Smith 
536aa482453SBarry Smith #if defined(PETSC_USE_LOG)
53777431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N);
5388965ea79SLois Curfman McInnes #endif
5398798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
540606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
5413501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
5423501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
5433501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
544622d7880SLois Curfman McInnes   if (mdn->factor) {
545606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
546606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
547606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
548606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
549622d7880SLois Curfman McInnes   }
550606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
551901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
552901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5533a40ed3dSBarry Smith   PetscFunctionReturn(0);
5548965ea79SLois Curfman McInnes }
55539ddd567SLois Curfman McInnes 
5564a2ae208SSatish Balay #undef __FUNCT__
5574a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5586849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5598965ea79SLois Curfman McInnes {
56039ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
561dfbe8321SBarry Smith   PetscErrorCode ierr;
5627056b6fcSBarry Smith 
5633a40ed3dSBarry Smith   PetscFunctionBegin;
56439ddd567SLois Curfman McInnes   if (mdn->size == 1) {
56539ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5668965ea79SLois Curfman McInnes   }
56729bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5683a40ed3dSBarry Smith   PetscFunctionReturn(0);
5698965ea79SLois Curfman McInnes }
5708965ea79SLois Curfman McInnes 
5714a2ae208SSatish Balay #undef __FUNCT__
5724a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
5736849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5748965ea79SLois Curfman McInnes {
57539ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
576dfbe8321SBarry Smith   PetscErrorCode    ierr;
57732dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
578b0a32e0cSBarry Smith   PetscViewerType   vtype;
57932077d6dSBarry Smith   PetscTruth        iascii,isdraw;
580b0a32e0cSBarry Smith   PetscViewer       sviewer;
581f3ef73ceSBarry Smith   PetscViewerFormat format;
5828965ea79SLois Curfman McInnes 
5833a40ed3dSBarry Smith   PetscFunctionBegin;
58432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
585fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
58632077d6dSBarry Smith   if (iascii) {
587b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
588b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
589456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
5904e220ebcSLois Curfman McInnes       MatInfo info;
591888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
59277431f27SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m,
59377431f27SBarry Smith                    (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr);
594b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5953501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5963a40ed3dSBarry Smith       PetscFunctionReturn(0);
597fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5983a40ed3dSBarry Smith       PetscFunctionReturn(0);
5998965ea79SLois Curfman McInnes     }
600f1af5d2fSBarry Smith   } else if (isdraw) {
601b0a32e0cSBarry Smith     PetscDraw       draw;
602f1af5d2fSBarry Smith     PetscTruth isnull;
603f1af5d2fSBarry Smith 
604b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
605b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
606f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
607f1af5d2fSBarry Smith   }
60877ed5343SBarry Smith 
6098965ea79SLois Curfman McInnes   if (size == 1) {
61039ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6113a40ed3dSBarry Smith   } else {
6128965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
6138965ea79SLois Curfman McInnes     Mat         A;
61413f74950SBarry Smith     PetscInt    M = mat->M,N = mat->N,m,row,i,nz;
615ba8c8a56SBarry Smith     PetscInt    *cols;
616ba8c8a56SBarry Smith     PetscScalar *vals;
6178965ea79SLois Curfman McInnes 
618f69a0ea3SMatthew Knepley     ierr = MatCreate(mat->comm,&A);CHKERRQ(ierr);
6198965ea79SLois Curfman McInnes     if (!rank) {
620f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
6213a40ed3dSBarry Smith     } else {
622f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
6238965ea79SLois Curfman McInnes     }
624be5d1d56SKris Buschelman     /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */
625878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
626878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
62752e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
6288965ea79SLois Curfman McInnes 
62939ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
63039ddd567SLois Curfman McInnes        but it's quick for now */
63151022da4SBarry Smith     A->insertmode = INSERT_VALUES;
632273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
6338965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
634ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
635ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
636ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
63739ddd567SLois Curfman McInnes       row++;
6388965ea79SLois Curfman McInnes     }
6398965ea79SLois Curfman McInnes 
6406d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6416d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
642b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
643b9b97703SBarry Smith     if (!rank) {
6446831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
6458965ea79SLois Curfman McInnes     }
646b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
647b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6488965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
6498965ea79SLois Curfman McInnes   }
6503a40ed3dSBarry Smith   PetscFunctionReturn(0);
6518965ea79SLois Curfman McInnes }
6528965ea79SLois Curfman McInnes 
6534a2ae208SSatish Balay #undef __FUNCT__
6544a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
655dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
6568965ea79SLois Curfman McInnes {
657dfbe8321SBarry Smith   PetscErrorCode ierr;
65832077d6dSBarry Smith   PetscTruth     iascii,isbinary,isdraw,issocket;
6598965ea79SLois Curfman McInnes 
660433994e6SBarry Smith   PetscFunctionBegin;
6610f5bd95cSBarry Smith 
66232077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
663fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
664b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
665fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6660f5bd95cSBarry Smith 
66732077d6dSBarry Smith   if (iascii || issocket || isdraw) {
668f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6690f5bd95cSBarry Smith   } else if (isbinary) {
6703a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6715cd90555SBarry Smith   } else {
672958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6738965ea79SLois Curfman McInnes   }
6743a40ed3dSBarry Smith   PetscFunctionReturn(0);
6758965ea79SLois Curfman McInnes }
6768965ea79SLois Curfman McInnes 
6774a2ae208SSatish Balay #undef __FUNCT__
6784a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
679dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6808965ea79SLois Curfman McInnes {
6813501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
6823501a2bdSLois Curfman McInnes   Mat            mdn = mat->A;
683dfbe8321SBarry Smith   PetscErrorCode ierr;
684329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
6858965ea79SLois Curfman McInnes 
6863a40ed3dSBarry Smith   PetscFunctionBegin;
687273d9f13SBarry Smith   info->rows_global    = (double)A->M;
688273d9f13SBarry Smith   info->columns_global = (double)A->N;
689273d9f13SBarry Smith   info->rows_local     = (double)A->m;
690273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6914e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6924e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6934e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6944e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6958965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6964e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6974e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6984e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6994e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7004e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7018965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
702d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
7034e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7044e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7054e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7064e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7074e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7088965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
709d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
7104e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7114e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7124e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7134e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7144e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7158965ea79SLois Curfman McInnes   }
7164e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
7174e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
7184e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
7193a40ed3dSBarry Smith   PetscFunctionReturn(0);
7208965ea79SLois Curfman McInnes }
7218965ea79SLois Curfman McInnes 
7224a2ae208SSatish Balay #undef __FUNCT__
7234a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
724dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op)
7258965ea79SLois Curfman McInnes {
72639ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
727dfbe8321SBarry Smith   PetscErrorCode ierr;
7288965ea79SLois Curfman McInnes 
7293a40ed3dSBarry Smith   PetscFunctionBegin;
73012c028f9SKris Buschelman   switch (op) {
73112c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
73212c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
73312c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
73412c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
73512c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
73612c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
737273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
73812c028f9SKris Buschelman     break;
73912c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
740273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
741273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
74212c028f9SKris Buschelman     break;
74312c028f9SKris Buschelman   case MAT_ROWS_SORTED:
74412c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
74512c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
74612c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
74763ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatSetOption_MPIDense:Option ignored\n"));CHKERRQ(ierr);
74812c028f9SKris Buschelman     break;
74912c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
750273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
751273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
75212c028f9SKris Buschelman     break;
75312c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
754273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
75512c028f9SKris Buschelman     break;
75612c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
75729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
75877e54ba9SKris Buschelman   case MAT_SYMMETRIC:
75977e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
7609a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
7619a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
7629a4540c5SBarry Smith   case MAT_HERMITIAN:
7639a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
7649a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
7659a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
76677e54ba9SKris Buschelman     break;
76712c028f9SKris Buschelman   default:
76829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
7693a40ed3dSBarry Smith   }
7703a40ed3dSBarry Smith   PetscFunctionReturn(0);
7718965ea79SLois Curfman McInnes }
7728965ea79SLois Curfman McInnes 
7738965ea79SLois Curfman McInnes 
7744a2ae208SSatish Balay #undef __FUNCT__
7754a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
776dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7775b2fa520SLois Curfman McInnes {
7785b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
7795b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
78087828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
781dfbe8321SBarry Smith   PetscErrorCode ierr;
78213f74950SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7835b2fa520SLois Curfman McInnes 
7845b2fa520SLois Curfman McInnes   PetscFunctionBegin;
78572d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7865b2fa520SLois Curfman McInnes   if (ll) {
78772d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
78829bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
7891ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
7905b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7915b2fa520SLois Curfman McInnes       x = l[i];
7925b2fa520SLois Curfman McInnes       v = mat->v + i;
7935b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7945b2fa520SLois Curfman McInnes     }
7951ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
796efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
7975b2fa520SLois Curfman McInnes   }
7985b2fa520SLois Curfman McInnes   if (rr) {
79972d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
80029bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
8015b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8025b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8031ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8045b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8055b2fa520SLois Curfman McInnes       x = r[i];
8065b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8075b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
8085b2fa520SLois Curfman McInnes     }
8091ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
810efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8115b2fa520SLois Curfman McInnes   }
8125b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8135b2fa520SLois Curfman McInnes }
8145b2fa520SLois Curfman McInnes 
8154a2ae208SSatish Balay #undef __FUNCT__
8164a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
817dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
818096963f5SLois Curfman McInnes {
8193501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8203501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
821dfbe8321SBarry Smith   PetscErrorCode ierr;
82213f74950SBarry Smith   PetscInt       i,j;
823329f5518SBarry Smith   PetscReal      sum = 0.0;
82487828ca2SBarry Smith   PetscScalar    *v = mat->v;
8253501a2bdSLois Curfman McInnes 
8263a40ed3dSBarry Smith   PetscFunctionBegin;
8273501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
828064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
8293501a2bdSLois Curfman McInnes   } else {
8303501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
831273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
832aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
833329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
8343501a2bdSLois Curfman McInnes #else
8353501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
8363501a2bdSLois Curfman McInnes #endif
8373501a2bdSLois Curfman McInnes       }
838064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
839064f8208SBarry Smith       *nrm = sqrt(*nrm);
840efee365bSSatish Balay       ierr = PetscLogFlops(2*mdn->A->n*mdn->A->m);CHKERRQ(ierr);
8413a40ed3dSBarry Smith     } else if (type == NORM_1) {
842329f5518SBarry Smith       PetscReal *tmp,*tmp2;
843b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
844273d9f13SBarry Smith       tmp2 = tmp + A->N;
845273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
846064f8208SBarry Smith       *nrm = 0.0;
8473501a2bdSLois Curfman McInnes       v = mat->v;
848273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
849273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
85067e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8513501a2bdSLois Curfman McInnes         }
8523501a2bdSLois Curfman McInnes       }
853d7d1e502SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
854273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
855064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8563501a2bdSLois Curfman McInnes       }
857606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
858efee365bSSatish Balay       ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
8593a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
860329f5518SBarry Smith       PetscReal ntemp;
8613501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
862064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
8633a40ed3dSBarry Smith     } else {
86429bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8653501a2bdSLois Curfman McInnes     }
8663501a2bdSLois Curfman McInnes   }
8673a40ed3dSBarry Smith   PetscFunctionReturn(0);
8683501a2bdSLois Curfman McInnes }
8693501a2bdSLois Curfman McInnes 
8704a2ae208SSatish Balay #undef __FUNCT__
8714a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
872dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout)
8733501a2bdSLois Curfman McInnes {
8743501a2bdSLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
8753501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
8763501a2bdSLois Curfman McInnes   Mat            B;
87713f74950SBarry Smith   PetscInt       M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8786849ba73SBarry Smith   PetscErrorCode ierr;
87913f74950SBarry Smith   PetscInt       j,i;
88087828ca2SBarry Smith   PetscScalar    *v;
8813501a2bdSLois Curfman McInnes 
8823a40ed3dSBarry Smith   PetscFunctionBegin;
8837c922b88SBarry Smith   if (!matout && M != N) {
88429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8857056b6fcSBarry Smith   }
886f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
887f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr);
888878740d9SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
889878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
8903501a2bdSLois Curfman McInnes 
891273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
89213f74950SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
8933501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8943501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8953501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8963501a2bdSLois Curfman McInnes     v   += m;
8973501a2bdSLois Curfman McInnes   }
898606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8996d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9006d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9017c922b88SBarry Smith   if (matout) {
9023501a2bdSLois Curfman McInnes     *matout = B;
9033501a2bdSLois Curfman McInnes   } else {
904273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
9053501a2bdSLois Curfman McInnes   }
9063a40ed3dSBarry Smith   PetscFunctionReturn(0);
907096963f5SLois Curfman McInnes }
908096963f5SLois Curfman McInnes 
909d9eff348SSatish Balay #include "petscblaslapack.h"
9104a2ae208SSatish Balay #undef __FUNCT__
9114a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
912f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha)
91344cd7ae7SLois Curfman McInnes {
91444cd7ae7SLois Curfman McInnes   Mat_MPIDense   *A = (Mat_MPIDense*)inA->data;
91544cd7ae7SLois Curfman McInnes   Mat_SeqDense   *a = (Mat_SeqDense*)A->A->data;
916f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
9174ce68768SBarry Smith   PetscBLASInt   one = 1,nz = (PetscBLASInt)inA->m*inA->N;
918efee365bSSatish Balay   PetscErrorCode ierr;
91944cd7ae7SLois Curfman McInnes 
9203a40ed3dSBarry Smith   PetscFunctionBegin;
921f4df32b1SMatthew Knepley   BLASscal_(&nz,&oalpha,a->v,&one);
922efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
9233a40ed3dSBarry Smith   PetscFunctionReturn(0);
92444cd7ae7SLois Curfman McInnes }
92544cd7ae7SLois Curfman McInnes 
9266849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
9278965ea79SLois Curfman McInnes 
9284a2ae208SSatish Balay #undef __FUNCT__
9294a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
930dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
931273d9f13SBarry Smith {
932dfbe8321SBarry Smith   PetscErrorCode ierr;
933273d9f13SBarry Smith 
934273d9f13SBarry Smith   PetscFunctionBegin;
935273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
936273d9f13SBarry Smith   PetscFunctionReturn(0);
937273d9f13SBarry Smith }
938273d9f13SBarry Smith 
9398965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
94009dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
94109dc0095SBarry Smith        MatGetRow_MPIDense,
94209dc0095SBarry Smith        MatRestoreRow_MPIDense,
94309dc0095SBarry Smith        MatMult_MPIDense,
94497304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
9457c922b88SBarry Smith        MatMultTranspose_MPIDense,
9467c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9478965ea79SLois Curfman McInnes        0,
94809dc0095SBarry Smith        0,
94909dc0095SBarry Smith        0,
95097304618SKris Buschelman /*10*/ 0,
95109dc0095SBarry Smith        0,
95209dc0095SBarry Smith        0,
95309dc0095SBarry Smith        0,
95409dc0095SBarry Smith        MatTranspose_MPIDense,
95597304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
95697304618SKris Buschelman        0,
95709dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9585b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
95909dc0095SBarry Smith        MatNorm_MPIDense,
96097304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
96109dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
96209dc0095SBarry Smith        0,
96309dc0095SBarry Smith        MatSetOption_MPIDense,
96409dc0095SBarry Smith        MatZeroEntries_MPIDense,
96597304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
96609dc0095SBarry Smith        0,
96709dc0095SBarry Smith        0,
96809dc0095SBarry Smith        0,
96909dc0095SBarry Smith        0,
97097304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
971273d9f13SBarry Smith        0,
97209dc0095SBarry Smith        0,
97309dc0095SBarry Smith        MatGetArray_MPIDense,
97409dc0095SBarry Smith        MatRestoreArray_MPIDense,
97597304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
97609dc0095SBarry Smith        0,
97709dc0095SBarry Smith        0,
97809dc0095SBarry Smith        0,
97909dc0095SBarry Smith        0,
98097304618SKris Buschelman /*40*/ 0,
9812ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
98209dc0095SBarry Smith        0,
98309dc0095SBarry Smith        MatGetValues_MPIDense,
98409dc0095SBarry Smith        0,
98597304618SKris Buschelman /*45*/ 0,
98609dc0095SBarry Smith        MatScale_MPIDense,
98709dc0095SBarry Smith        0,
98809dc0095SBarry Smith        0,
98909dc0095SBarry Smith        0,
990521d7252SBarry Smith /*50*/ 0,
99109dc0095SBarry Smith        0,
99209dc0095SBarry Smith        0,
99309dc0095SBarry Smith        0,
99409dc0095SBarry Smith        0,
99597304618SKris Buschelman /*55*/ 0,
99609dc0095SBarry Smith        0,
99709dc0095SBarry Smith        0,
99809dc0095SBarry Smith        0,
99909dc0095SBarry Smith        0,
100097304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
1001b9b97703SBarry Smith        MatDestroy_MPIDense,
1002b9b97703SBarry Smith        MatView_MPIDense,
100397304618SKris Buschelman        MatGetPetscMaps_Petsc,
100497304618SKris Buschelman        0,
100597304618SKris Buschelman /*65*/ 0,
100697304618SKris Buschelman        0,
100797304618SKris Buschelman        0,
100897304618SKris Buschelman        0,
100997304618SKris Buschelman        0,
101097304618SKris Buschelman /*70*/ 0,
101197304618SKris Buschelman        0,
101297304618SKris Buschelman        0,
101397304618SKris Buschelman        0,
101497304618SKris Buschelman        0,
101597304618SKris Buschelman /*75*/ 0,
101697304618SKris Buschelman        0,
101797304618SKris Buschelman        0,
101897304618SKris Buschelman        0,
101997304618SKris Buschelman        0,
102097304618SKris Buschelman /*80*/ 0,
102197304618SKris Buschelman        0,
102297304618SKris Buschelman        0,
102397304618SKris Buschelman        0,
1024865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense,
1025865e5f61SKris Buschelman        0,
1026865e5f61SKris Buschelman        0,
1027865e5f61SKris Buschelman        0,
1028865e5f61SKris Buschelman        0,
1029865e5f61SKris Buschelman        0,
1030865e5f61SKris Buschelman /*90*/ 0,
1031865e5f61SKris Buschelman        0,
1032865e5f61SKris Buschelman        0,
1033865e5f61SKris Buschelman        0,
1034865e5f61SKris Buschelman        0,
1035865e5f61SKris Buschelman /*95*/ 0,
1036865e5f61SKris Buschelman        0,
1037865e5f61SKris Buschelman        0,
1038865e5f61SKris Buschelman        0};
10398965ea79SLois Curfman McInnes 
1040273d9f13SBarry Smith EXTERN_C_BEGIN
10414a2ae208SSatish Balay #undef __FUNCT__
1042a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
1043be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1044a23d5eceSKris Buschelman {
1045a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1046dfbe8321SBarry Smith   PetscErrorCode ierr;
1047a23d5eceSKris Buschelman 
1048a23d5eceSKris Buschelman   PetscFunctionBegin;
1049a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1050a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1051a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1052a23d5eceSKris Buschelman 
1053a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
1054f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1055f69a0ea3SMatthew Knepley   ierr = MatSetSizes(a->A,mat->m,mat->N,mat->m,mat->N);CHKERRQ(ierr);
10565c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
10575c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
105852e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1059a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1060a23d5eceSKris Buschelman }
1061a23d5eceSKris Buschelman EXTERN_C_END
1062a23d5eceSKris Buschelman 
10630bad9183SKris Buschelman /*MC
1064fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
10650bad9183SKris Buschelman 
10660bad9183SKris Buschelman    Options Database Keys:
10670bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
10680bad9183SKris Buschelman 
10690bad9183SKris Buschelman   Level: beginner
10700bad9183SKris Buschelman 
10710bad9183SKris Buschelman .seealso: MatCreateMPIDense
10720bad9183SKris Buschelman M*/
10730bad9183SKris Buschelman 
1074a23d5eceSKris Buschelman EXTERN_C_BEGIN
1075a23d5eceSKris Buschelman #undef __FUNCT__
10764a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1077be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat)
1078273d9f13SBarry Smith {
1079273d9f13SBarry Smith   Mat_MPIDense   *a;
1080dfbe8321SBarry Smith   PetscErrorCode ierr;
108113f74950SBarry Smith   PetscInt       i;
1082273d9f13SBarry Smith 
1083273d9f13SBarry Smith   PetscFunctionBegin;
1084b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1085b0a32e0cSBarry Smith   mat->data         = (void*)a;
1086273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1087273d9f13SBarry Smith   mat->factor       = 0;
1088273d9f13SBarry Smith   mat->mapping      = 0;
1089273d9f13SBarry Smith 
1090273d9f13SBarry Smith   a->factor       = 0;
1091273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1092273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
1093273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
1094273d9f13SBarry Smith 
1095273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
1096273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
1097273d9f13SBarry Smith   a->nvec = mat->n;
1098273d9f13SBarry Smith 
1099273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1100273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
11018a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
11028a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
1103273d9f13SBarry Smith 
1104273d9f13SBarry Smith   /* build local table of row and column ownerships */
110513f74950SBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr);
1106273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
110752e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr);
110813f74950SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr);
1109273d9f13SBarry Smith   a->rowners[0] = 0;
1110273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1111273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
1112273d9f13SBarry Smith   }
1113273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
1114273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
111513f74950SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr);
1116273d9f13SBarry Smith   a->cowners[0] = 0;
1117273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1118273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1119273d9f13SBarry Smith   }
1120273d9f13SBarry Smith 
1121273d9f13SBarry Smith   /* build cache for off array entries formed */
1122273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1123273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1124273d9f13SBarry Smith 
1125273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1126273d9f13SBarry Smith   a->lvec        = 0;
1127273d9f13SBarry Smith   a->Mvctx       = 0;
1128273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1129273d9f13SBarry Smith 
1130273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1131273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1132273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1133a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1134a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1135a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1136273d9f13SBarry Smith   PetscFunctionReturn(0);
1137273d9f13SBarry Smith }
1138273d9f13SBarry Smith EXTERN_C_END
1139273d9f13SBarry Smith 
1140209238afSKris Buschelman /*MC
1141002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1142209238afSKris Buschelman 
1143209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1144209238afSKris Buschelman    and MATMPIDENSE otherwise.
1145209238afSKris Buschelman 
1146209238afSKris Buschelman    Options Database Keys:
1147209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1148209238afSKris Buschelman 
1149209238afSKris Buschelman   Level: beginner
1150209238afSKris Buschelman 
1151209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1152209238afSKris Buschelman M*/
1153209238afSKris Buschelman 
1154209238afSKris Buschelman EXTERN_C_BEGIN
1155209238afSKris Buschelman #undef __FUNCT__
1156209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1157be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A)
1158dfbe8321SBarry Smith {
11596849ba73SBarry Smith   PetscErrorCode ierr;
116013f74950SBarry Smith   PetscMPIInt    size;
1161209238afSKris Buschelman 
1162209238afSKris Buschelman   PetscFunctionBegin;
1163209238afSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr);
1164209238afSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1165209238afSKris Buschelman   if (size == 1) {
1166209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1167209238afSKris Buschelman   } else {
1168209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1169209238afSKris Buschelman   }
1170209238afSKris Buschelman   PetscFunctionReturn(0);
1171209238afSKris Buschelman }
1172209238afSKris Buschelman EXTERN_C_END
1173209238afSKris Buschelman 
11744a2ae208SSatish Balay #undef __FUNCT__
11754a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1176273d9f13SBarry Smith /*@C
1177273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1178273d9f13SBarry Smith 
1179273d9f13SBarry Smith    Not collective
1180273d9f13SBarry Smith 
1181273d9f13SBarry Smith    Input Parameters:
1182273d9f13SBarry Smith .  A - the matrix
1183273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1184273d9f13SBarry Smith    to control all matrix memory allocation.
1185273d9f13SBarry Smith 
1186273d9f13SBarry Smith    Notes:
1187273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1188273d9f13SBarry Smith    storage by columns.
1189273d9f13SBarry Smith 
1190273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1191273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1192273d9f13SBarry Smith    set data=PETSC_NULL.
1193273d9f13SBarry Smith 
1194273d9f13SBarry Smith    Level: intermediate
1195273d9f13SBarry Smith 
1196273d9f13SBarry Smith .keywords: matrix,dense, parallel
1197273d9f13SBarry Smith 
1198273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1199273d9f13SBarry Smith @*/
1200be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1201273d9f13SBarry Smith {
1202dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar *);
1203273d9f13SBarry Smith 
1204273d9f13SBarry Smith   PetscFunctionBegin;
1205565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1206a23d5eceSKris Buschelman   if (f) {
1207a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1208a23d5eceSKris Buschelman   }
1209273d9f13SBarry Smith   PetscFunctionReturn(0);
1210273d9f13SBarry Smith }
1211273d9f13SBarry Smith 
12124a2ae208SSatish Balay #undef __FUNCT__
12134a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
12148965ea79SLois Curfman McInnes /*@C
121539ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
12168965ea79SLois Curfman McInnes 
1217db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1218db81eaa0SLois Curfman McInnes 
12198965ea79SLois Curfman McInnes    Input Parameters:
1220db81eaa0SLois Curfman McInnes +  comm - MPI communicator
12218965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1222db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
12238965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1224db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
12257f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1226dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
12278965ea79SLois Curfman McInnes 
12288965ea79SLois Curfman McInnes    Output Parameter:
1229477f1c0bSLois Curfman McInnes .  A - the matrix
12308965ea79SLois Curfman McInnes 
1231b259b22eSLois Curfman McInnes    Notes:
123239ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
123339ddd567SLois Curfman McInnes    storage by columns.
12348965ea79SLois Curfman McInnes 
123518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
123618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
12377f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
123818f449edSLois Curfman McInnes 
12398965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
12408965ea79SLois Curfman McInnes    (possibly both).
12418965ea79SLois Curfman McInnes 
1242027ccd11SLois Curfman McInnes    Level: intermediate
1243027ccd11SLois Curfman McInnes 
124439ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
12458965ea79SLois Curfman McInnes 
124639ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
12478965ea79SLois Curfman McInnes @*/
1248be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
12498965ea79SLois Curfman McInnes {
12506849ba73SBarry Smith   PetscErrorCode ierr;
125113f74950SBarry Smith   PetscMPIInt    size;
12528965ea79SLois Curfman McInnes 
12533a40ed3dSBarry Smith   PetscFunctionBegin;
1254f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1255f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1256273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1257273d9f13SBarry Smith   if (size > 1) {
1258273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1259273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1260273d9f13SBarry Smith   } else {
1261273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1262273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
12638c469469SLois Curfman McInnes   }
12643a40ed3dSBarry Smith   PetscFunctionReturn(0);
12658965ea79SLois Curfman McInnes }
12668965ea79SLois Curfman McInnes 
12674a2ae208SSatish Balay #undef __FUNCT__
12684a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
12696849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12708965ea79SLois Curfman McInnes {
12718965ea79SLois Curfman McInnes   Mat            mat;
12723501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1273dfbe8321SBarry Smith   PetscErrorCode ierr;
12748965ea79SLois Curfman McInnes 
12753a40ed3dSBarry Smith   PetscFunctionBegin;
12768965ea79SLois Curfman McInnes   *newmat       = 0;
1277f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&mat);CHKERRQ(ierr);
1278f69a0ea3SMatthew Knepley   ierr = MatSetSizes(mat,A->m,A->n,A->M,A->N);CHKERRQ(ierr);
1279be5d1d56SKris Buschelman   ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr);
1280b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1281b0a32e0cSBarry Smith   mat->data         = (void*)a;
1282549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
12833501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1284c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1285273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
12868965ea79SLois Curfman McInnes 
12878965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
12888965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
12898965ea79SLois Curfman McInnes   a->size         = oldmat->size;
12908965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1291e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1292b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
12933782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
129413f74950SBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr);
129552e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr);
129613f74950SBarry Smith   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
12978798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
12988965ea79SLois Curfman McInnes 
1299329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
13005609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
130152e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
13028965ea79SLois Curfman McInnes   *newmat = mat;
13033a40ed3dSBarry Smith   PetscFunctionReturn(0);
13048965ea79SLois Curfman McInnes }
13058965ea79SLois Curfman McInnes 
1306e090d566SSatish Balay #include "petscsys.h"
13078965ea79SLois Curfman McInnes 
13084a2ae208SSatish Balay #undef __FUNCT__
13094a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
1310f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat)
131190ace30eSBarry Smith {
13126849ba73SBarry Smith   PetscErrorCode ierr;
131332dcc486SBarry Smith   PetscMPIInt    rank,size;
131413f74950SBarry Smith   PetscInt       *rowners,i,m,nz,j;
131587828ca2SBarry Smith   PetscScalar    *array,*vals,*vals_ptr;
131690ace30eSBarry Smith   MPI_Status     status;
131790ace30eSBarry Smith 
13183a40ed3dSBarry Smith   PetscFunctionBegin;
1319d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1320d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
132190ace30eSBarry Smith 
132290ace30eSBarry Smith   /* determine ownership of all rows */
132390ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
132413f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
132513f74950SBarry Smith   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
132690ace30eSBarry Smith   rowners[0] = 0;
132790ace30eSBarry Smith   for (i=2; i<=size; i++) {
132890ace30eSBarry Smith     rowners[i] += rowners[i-1];
132990ace30eSBarry Smith   }
133090ace30eSBarry Smith 
1331f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1332f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1333be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1334878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
133590ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
133690ace30eSBarry Smith 
133790ace30eSBarry Smith   if (!rank) {
133887828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
133990ace30eSBarry Smith 
134090ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
13410752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
134290ace30eSBarry Smith 
134390ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
134490ace30eSBarry Smith     vals_ptr = vals;
134590ace30eSBarry Smith     for (i=0; i<m; i++) {
134690ace30eSBarry Smith       for (j=0; j<N; j++) {
134790ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
134890ace30eSBarry Smith       }
134990ace30eSBarry Smith     }
135090ace30eSBarry Smith 
135190ace30eSBarry Smith     /* read in other processors and ship out */
135290ace30eSBarry Smith     for (i=1; i<size; i++) {
135390ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
13540752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1355ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
135690ace30eSBarry Smith     }
13573a40ed3dSBarry Smith   } else {
135890ace30eSBarry Smith     /* receive numeric values */
135987828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
136090ace30eSBarry Smith 
136190ace30eSBarry Smith     /* receive message of values*/
1362ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
136390ace30eSBarry Smith 
136490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
136590ace30eSBarry Smith     vals_ptr = vals;
136690ace30eSBarry Smith     for (i=0; i<m; i++) {
136790ace30eSBarry Smith       for (j=0; j<N; j++) {
136890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
136990ace30eSBarry Smith       }
137090ace30eSBarry Smith     }
137190ace30eSBarry Smith   }
1372606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1373606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
13746d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13756d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13763a40ed3dSBarry Smith   PetscFunctionReturn(0);
137790ace30eSBarry Smith }
137890ace30eSBarry Smith 
13794a2ae208SSatish Balay #undef __FUNCT__
13804a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
1381f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat)
13828965ea79SLois Curfman McInnes {
13838965ea79SLois Curfman McInnes   Mat            A;
138487828ca2SBarry Smith   PetscScalar    *vals,*svals;
138519bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
13868965ea79SLois Curfman McInnes   MPI_Status     status;
138713f74950SBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
138813f74950SBarry Smith   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
138913f74950SBarry Smith   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
139013f74950SBarry Smith   PetscInt       i,nz,j,rstart,rend;
139113f74950SBarry Smith   int            fd;
13926849ba73SBarry Smith   PetscErrorCode ierr;
13938965ea79SLois Curfman McInnes 
13943a40ed3dSBarry Smith   PetscFunctionBegin;
1395d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1396d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13978965ea79SLois Curfman McInnes   if (!rank) {
1398b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
13990752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1400552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
14018965ea79SLois Curfman McInnes   }
14028965ea79SLois Curfman McInnes 
140313f74950SBarry Smith   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
140490ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
140590ace30eSBarry Smith 
140690ace30eSBarry Smith   /*
140790ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
140890ace30eSBarry Smith   */
140990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1410be5d1d56SKris Buschelman     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr);
14113a40ed3dSBarry Smith     PetscFunctionReturn(0);
141290ace30eSBarry Smith   }
141390ace30eSBarry Smith 
14148965ea79SLois Curfman McInnes   /* determine ownership of all rows */
14158965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
141613f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
1417ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
14188965ea79SLois Curfman McInnes   rowners[0] = 0;
14198965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
14208965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
14218965ea79SLois Curfman McInnes   }
14228965ea79SLois Curfman McInnes   rstart = rowners[rank];
14238965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
14248965ea79SLois Curfman McInnes 
14258965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
142613f74950SBarry Smith   ierr    = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr);
14278965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
14288965ea79SLois Curfman McInnes   if (!rank) {
142913f74950SBarry Smith     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
14300752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
143113f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
14328965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
14331466f1e1SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
1434606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1435ca161407SBarry Smith   } else {
14361466f1e1SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
14378965ea79SLois Curfman McInnes   }
14388965ea79SLois Curfman McInnes 
14398965ea79SLois Curfman McInnes   if (!rank) {
14408965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
144113f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
144213f74950SBarry Smith     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
14438965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14448965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
14458965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
14468965ea79SLois Curfman McInnes       }
14478965ea79SLois Curfman McInnes     }
1448606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
14498965ea79SLois Curfman McInnes 
14508965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
14518965ea79SLois Curfman McInnes     maxnz = 0;
14528965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14530452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
14548965ea79SLois Curfman McInnes     }
145513f74950SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
14568965ea79SLois Curfman McInnes 
14578965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
14588965ea79SLois Curfman McInnes     nz = procsnz[0];
145913f74950SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
14600752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
14618965ea79SLois Curfman McInnes 
14628965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
14638965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14648965ea79SLois Curfman McInnes       nz   = procsnz[i];
14650752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
146613f74950SBarry Smith       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
14678965ea79SLois Curfman McInnes     }
1468606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
14693a40ed3dSBarry Smith   } else {
14708965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
14718965ea79SLois Curfman McInnes     nz = 0;
14728965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14738965ea79SLois Curfman McInnes       nz += ourlens[i];
14748965ea79SLois Curfman McInnes     }
147513f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
14768965ea79SLois Curfman McInnes 
14778965ea79SLois Curfman McInnes     /* receive message of column indices*/
147813f74950SBarry Smith     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
147913f74950SBarry Smith     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
148029bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14818965ea79SLois Curfman McInnes   }
14828965ea79SLois Curfman McInnes 
14838965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
148413f74950SBarry Smith   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
14858965ea79SLois Curfman McInnes   jj = 0;
14868965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14878965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
14888965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
14898965ea79SLois Curfman McInnes       jj++;
14908965ea79SLois Curfman McInnes     }
14918965ea79SLois Curfman McInnes   }
14928965ea79SLois Curfman McInnes 
14938965ea79SLois Curfman McInnes   /* create our matrix */
14948965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14958965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
14968965ea79SLois Curfman McInnes   }
1497f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1498f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1499be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1500878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
15018965ea79SLois Curfman McInnes   A = *newmat;
15028965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
15038965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
15048965ea79SLois Curfman McInnes   }
15058965ea79SLois Curfman McInnes 
15068965ea79SLois Curfman McInnes   if (!rank) {
150787828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
15088965ea79SLois Curfman McInnes 
15098965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
15108965ea79SLois Curfman McInnes     nz = procsnz[0];
15110752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
15128965ea79SLois Curfman McInnes 
15138965ea79SLois Curfman McInnes     /* insert into matrix */
15148965ea79SLois Curfman McInnes     jj      = rstart;
15158965ea79SLois Curfman McInnes     smycols = mycols;
15168965ea79SLois Curfman McInnes     svals   = vals;
15178965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15188965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15198965ea79SLois Curfman McInnes       smycols += ourlens[i];
15208965ea79SLois Curfman McInnes       svals   += ourlens[i];
15218965ea79SLois Curfman McInnes       jj++;
15228965ea79SLois Curfman McInnes     }
15238965ea79SLois Curfman McInnes 
15248965ea79SLois Curfman McInnes     /* read in other processors and ship out */
15258965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
15268965ea79SLois Curfman McInnes       nz   = procsnz[i];
15270752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1528ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
15298965ea79SLois Curfman McInnes     }
1530606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
15313a40ed3dSBarry Smith   } else {
15328965ea79SLois Curfman McInnes     /* receive numeric values */
153384d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
15348965ea79SLois Curfman McInnes 
15358965ea79SLois Curfman McInnes     /* receive message of values*/
1536ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1537ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
153829bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
15398965ea79SLois Curfman McInnes 
15408965ea79SLois Curfman McInnes     /* insert into matrix */
15418965ea79SLois Curfman McInnes     jj      = rstart;
15428965ea79SLois Curfman McInnes     smycols = mycols;
15438965ea79SLois Curfman McInnes     svals   = vals;
15448965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15458965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15468965ea79SLois Curfman McInnes       smycols += ourlens[i];
15478965ea79SLois Curfman McInnes       svals   += ourlens[i];
15488965ea79SLois Curfman McInnes       jj++;
15498965ea79SLois Curfman McInnes     }
15508965ea79SLois Curfman McInnes   }
1551606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1552606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1553606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1554606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
15558965ea79SLois Curfman McInnes 
15566d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15576d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15583a40ed3dSBarry Smith   PetscFunctionReturn(0);
15598965ea79SLois Curfman McInnes }
156090ace30eSBarry Smith 
156190ace30eSBarry Smith 
156290ace30eSBarry Smith 
156390ace30eSBarry Smith 
156490ace30eSBarry Smith 
1565