xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 175be7b47374493d4cbd1171110ef131e9bbc2ec)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3ed3cc1f0SBarry Smith /*
4ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
5ed3cc1f0SBarry Smith */
6ed3cc1f0SBarry Smith 
789ac3759SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"    /*I   "petscmat.h"  I*/
88965ea79SLois Curfman McInnes 
9ba8c8a56SBarry Smith #undef __FUNCT__
100ad19415SHong Zhang #define __FUNCT__ "MatLUFactorSymbolic_MPIDense"
110ad19415SHong Zhang PetscErrorCode MatLUFactorSymbolic_MPIDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
120ad19415SHong Zhang {
130ad19415SHong Zhang   PetscErrorCode ierr;
140ad19415SHong Zhang 
150ad19415SHong Zhang   PetscFunctionBegin;
166017fc9fSHong Zhang   ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
170ad19415SHong Zhang   PetscFunctionReturn(0);
180ad19415SHong Zhang }
190ad19415SHong Zhang 
200ad19415SHong Zhang #undef __FUNCT__
210ad19415SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_MPIDense"
220ad19415SHong Zhang PetscErrorCode MatLUFactorNumeric_MPIDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
230ad19415SHong Zhang {
240ad19415SHong Zhang   PetscFunctionBegin;
250ad19415SHong Zhang   SETERRQ(PETSC_ERR_SUP,"No support of numerical factorization for mpidense matrix type");
260ad19415SHong Zhang   PetscFunctionReturn(0);
270ad19415SHong Zhang }
280ad19415SHong Zhang 
290ad19415SHong Zhang #undef __FUNCT__
300ad19415SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_MPIDense"
316017fc9fSHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_MPIDense(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
320ad19415SHong Zhang {
330ad19415SHong Zhang   PetscErrorCode ierr;
340ad19415SHong Zhang 
350ad19415SHong Zhang   PetscFunctionBegin;
366017fc9fSHong Zhang   ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
370ad19415SHong Zhang   PetscFunctionReturn(0);
380ad19415SHong Zhang }
390ad19415SHong Zhang 
400ad19415SHong Zhang #undef __FUNCT__
410ad19415SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_MPIDense"
420ad19415SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_MPIDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
430ad19415SHong Zhang {
440ad19415SHong Zhang   PetscFunctionBegin;
450ad19415SHong Zhang   SETERRQ(PETSC_ERR_SUP,"no support of numerical factorization for mpidense matrix type");
460ad19415SHong Zhang   PetscFunctionReturn(0);
470ad19415SHong Zhang }
480ad19415SHong Zhang 
490ad19415SHong Zhang #undef __FUNCT__
50ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
51ab92ecdeSBarry Smith /*@
52ab92ecdeSBarry Smith 
53ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
54ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
55ab92ecdeSBarry Smith 
56ab92ecdeSBarry Smith     Input Parameter:
57ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
58ab92ecdeSBarry Smith 
59ab92ecdeSBarry Smith     Output Parameter:
60ab92ecdeSBarry Smith .      B - the inner matrix
61ab92ecdeSBarry Smith 
628e6c10adSSatish Balay     Level: intermediate
638e6c10adSSatish Balay 
64ab92ecdeSBarry Smith @*/
65ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
66ab92ecdeSBarry Smith {
67ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
68ab92ecdeSBarry Smith   PetscErrorCode ierr;
69ab92ecdeSBarry Smith   PetscTruth     flg;
70ab92ecdeSBarry Smith   PetscMPIInt    size;
71ab92ecdeSBarry Smith 
72ab92ecdeSBarry Smith   PetscFunctionBegin;
73ab92ecdeSBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
74ab92ecdeSBarry Smith   if (!flg) {  /* this check sucks! */
75ab92ecdeSBarry Smith     ierr = PetscTypeCompare((PetscObject)A,MATDENSE,&flg);CHKERRQ(ierr);
76ab92ecdeSBarry Smith     if (flg) {
77ab92ecdeSBarry Smith       ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
78ab92ecdeSBarry Smith       if (size == 1) flg = PETSC_FALSE;
79ab92ecdeSBarry Smith     }
80ab92ecdeSBarry Smith   }
81ab92ecdeSBarry Smith   if (flg) {
82ab92ecdeSBarry Smith     *B = mat->A;
83ab92ecdeSBarry Smith   } else {
84ab92ecdeSBarry Smith     *B = A;
85ab92ecdeSBarry Smith   }
86ab92ecdeSBarry Smith   PetscFunctionReturn(0);
87ab92ecdeSBarry Smith }
88ab92ecdeSBarry Smith 
89ab92ecdeSBarry Smith #undef __FUNCT__
90ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
91ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
92ba8c8a56SBarry Smith {
93ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
94ba8c8a56SBarry Smith   PetscErrorCode ierr;
95ba8c8a56SBarry Smith   PetscInt       lrow,rstart = mat->rstart,rend = mat->rend;
96ba8c8a56SBarry Smith 
97ba8c8a56SBarry Smith   PetscFunctionBegin;
98ba8c8a56SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
99ba8c8a56SBarry Smith   lrow = row - rstart;
100ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
101ba8c8a56SBarry Smith   PetscFunctionReturn(0);
102ba8c8a56SBarry Smith }
103ba8c8a56SBarry Smith 
104ba8c8a56SBarry Smith #undef __FUNCT__
105ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
106ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
107ba8c8a56SBarry Smith {
108ba8c8a56SBarry Smith   PetscErrorCode ierr;
109ba8c8a56SBarry Smith 
110ba8c8a56SBarry Smith   PetscFunctionBegin;
111ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
112ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
113ba8c8a56SBarry Smith   PetscFunctionReturn(0);
114ba8c8a56SBarry Smith }
115ba8c8a56SBarry Smith 
1160de54da6SSatish Balay EXTERN_C_BEGIN
1174a2ae208SSatish Balay #undef __FUNCT__
1184a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
119be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
1200de54da6SSatish Balay {
1210de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
1226849ba73SBarry Smith   PetscErrorCode ierr;
12313f74950SBarry Smith   PetscInt       m = A->m,rstart = mdn->rstart;
12487828ca2SBarry Smith   PetscScalar    *array;
1250de54da6SSatish Balay   MPI_Comm       comm;
1260de54da6SSatish Balay 
1270de54da6SSatish Balay   PetscFunctionBegin;
128273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
1290de54da6SSatish Balay 
1300de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
1310de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
1320de54da6SSatish Balay 
1330de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
1340de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
135f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,B);CHKERRQ(ierr);
136f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr);
1375c5985e7SKris Buschelman   ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr);
1385c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
1390de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
1400de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1410de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1420de54da6SSatish Balay 
1430de54da6SSatish Balay   *iscopy = PETSC_TRUE;
1440de54da6SSatish Balay   PetscFunctionReturn(0);
1450de54da6SSatish Balay }
1460de54da6SSatish Balay EXTERN_C_END
1470de54da6SSatish Balay 
1484a2ae208SSatish Balay #undef __FUNCT__
1494a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
15013f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1518965ea79SLois Curfman McInnes {
15239b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
153dfbe8321SBarry Smith   PetscErrorCode ierr;
15413f74950SBarry Smith   PetscInt       i,j,rstart = A->rstart,rend = A->rend,row;
155273d9f13SBarry Smith   PetscTruth     roworiented = A->roworiented;
1568965ea79SLois Curfman McInnes 
1573a40ed3dSBarry Smith   PetscFunctionBegin;
1588965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1595ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
160273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1618965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1628965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
16339b7565bSBarry Smith       if (roworiented) {
16439b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1653a40ed3dSBarry Smith       } else {
1668965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1675ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
168273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
16939b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
17039b7565bSBarry Smith         }
1718965ea79SLois Curfman McInnes       }
1723a40ed3dSBarry Smith     } else {
1733782ba37SSatish Balay       if (!A->donotstash) {
17439b7565bSBarry Smith         if (roworiented) {
1758798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
176d36fbae8SSatish Balay         } else {
1778798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
17839b7565bSBarry Smith         }
179b49de8d1SLois Curfman McInnes       }
180b49de8d1SLois Curfman McInnes     }
1813782ba37SSatish Balay   }
1823a40ed3dSBarry Smith   PetscFunctionReturn(0);
183b49de8d1SLois Curfman McInnes }
184b49de8d1SLois Curfman McInnes 
1854a2ae208SSatish Balay #undef __FUNCT__
1864a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
18713f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
188b49de8d1SLois Curfman McInnes {
189b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
190dfbe8321SBarry Smith   PetscErrorCode ierr;
19113f74950SBarry Smith   PetscInt       i,j,rstart = mdn->rstart,rend = mdn->rend,row;
192b49de8d1SLois Curfman McInnes 
1933a40ed3dSBarry Smith   PetscFunctionBegin;
194b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
19529bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
196273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
197b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
198b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
199b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
20029bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
201273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
20229bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
203a8c6a408SBarry Smith         }
204b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
205b49de8d1SLois Curfman McInnes       }
206a8c6a408SBarry Smith     } else {
20729bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
2088965ea79SLois Curfman McInnes     }
2098965ea79SLois Curfman McInnes   }
2103a40ed3dSBarry Smith   PetscFunctionReturn(0);
2118965ea79SLois Curfman McInnes }
2128965ea79SLois Curfman McInnes 
2134a2ae208SSatish Balay #undef __FUNCT__
2144a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
215dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
216ff14e315SSatish Balay {
217ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
218dfbe8321SBarry Smith   PetscErrorCode ierr;
219ff14e315SSatish Balay 
2203a40ed3dSBarry Smith   PetscFunctionBegin;
221ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
2223a40ed3dSBarry Smith   PetscFunctionReturn(0);
223ff14e315SSatish Balay }
224ff14e315SSatish Balay 
2254a2ae208SSatish Balay #undef __FUNCT__
2264a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
22713f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
228ca3fa75bSLois Curfman McInnes {
229ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data,*newmatd;
230ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
2316849ba73SBarry Smith   PetscErrorCode ierr;
23213f74950SBarry Smith   PetscInt       i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
23387828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
234ca3fa75bSLois Curfman McInnes   Mat            newmat;
235ca3fa75bSLois Curfman McInnes 
236ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
237ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
238ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
239b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
240b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
241ca3fa75bSLois Curfman McInnes 
242ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
2437eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
2447eba5e9cSLois Curfman McInnes      original matrix! */
245ca3fa75bSLois Curfman McInnes 
246ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
2477eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
248ca3fa75bSLois Curfman McInnes 
249ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
250ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
25129bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2527eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
253ca3fa75bSLois Curfman McInnes     newmat = *B;
254ca3fa75bSLois Curfman McInnes   } else {
255ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
256f69a0ea3SMatthew Knepley     ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr);
257f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr);
258878740d9SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
259878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
260ca3fa75bSLois Curfman McInnes   }
261ca3fa75bSLois Curfman McInnes 
262ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
263ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
264ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
265ca3fa75bSLois Curfman McInnes 
266ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
267ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
268ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2697eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
270ca3fa75bSLois Curfman McInnes     }
271ca3fa75bSLois Curfman McInnes   }
272ca3fa75bSLois Curfman McInnes 
273ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
274ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
275ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
276ca3fa75bSLois Curfman McInnes 
277ca3fa75bSLois Curfman McInnes   /* Free work space */
278ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
279ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
280ca3fa75bSLois Curfman McInnes   *B = newmat;
281ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
282ca3fa75bSLois Curfman McInnes }
283ca3fa75bSLois Curfman McInnes 
2844a2ae208SSatish Balay #undef __FUNCT__
2854a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
286dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
287ff14e315SSatish Balay {
2883a40ed3dSBarry Smith   PetscFunctionBegin;
2893a40ed3dSBarry Smith   PetscFunctionReturn(0);
290ff14e315SSatish Balay }
291ff14e315SSatish Balay 
2924a2ae208SSatish Balay #undef __FUNCT__
2934a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
294dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2958965ea79SLois Curfman McInnes {
29639ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
2978965ea79SLois Curfman McInnes   MPI_Comm       comm = mat->comm;
298dfbe8321SBarry Smith   PetscErrorCode ierr;
29913f74950SBarry Smith   PetscInt       nstash,reallocs;
3008965ea79SLois Curfman McInnes   InsertMode     addv;
3018965ea79SLois Curfman McInnes 
3023a40ed3dSBarry Smith   PetscFunctionBegin;
3038965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
304ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
3057056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
30629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
3078965ea79SLois Curfman McInnes   }
308e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
3098965ea79SLois Curfman McInnes 
3108798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
3118798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
31263ba0a88SBarry Smith   ierr = PetscLogInfo((mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs));CHKERRQ(ierr);
3133a40ed3dSBarry Smith   PetscFunctionReturn(0);
3148965ea79SLois Curfman McInnes }
3158965ea79SLois Curfman McInnes 
3164a2ae208SSatish Balay #undef __FUNCT__
3174a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
318dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
3198965ea79SLois Curfman McInnes {
32039ddd567SLois Curfman McInnes   Mat_MPIDense    *mdn=(Mat_MPIDense*)mat->data;
3216849ba73SBarry Smith   PetscErrorCode  ierr;
32213f74950SBarry Smith   PetscInt        i,*row,*col,flg,j,rstart,ncols;
32313f74950SBarry Smith   PetscMPIInt     n;
32487828ca2SBarry Smith   PetscScalar     *val;
325e0fa3b82SLois Curfman McInnes   InsertMode      addv=mat->insertmode;
3268965ea79SLois Curfman McInnes 
3273a40ed3dSBarry Smith   PetscFunctionBegin;
3288965ea79SLois Curfman McInnes   /*  wait on receives */
3297ef1d9bdSSatish Balay   while (1) {
3308798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
3317ef1d9bdSSatish Balay     if (!flg) break;
3328965ea79SLois Curfman McInnes 
3337ef1d9bdSSatish Balay     for (i=0; i<n;) {
3347ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
3357ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
3367ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
3377ef1d9bdSSatish Balay       else       ncols = n-i;
3387ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
3397ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
3407ef1d9bdSSatish Balay       i = j;
3418965ea79SLois Curfman McInnes     }
3427ef1d9bdSSatish Balay   }
3438798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
3448965ea79SLois Curfman McInnes 
34539ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
34639ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3478965ea79SLois Curfman McInnes 
3486d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
34939ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3508965ea79SLois Curfman McInnes   }
3513a40ed3dSBarry Smith   PetscFunctionReturn(0);
3528965ea79SLois Curfman McInnes }
3538965ea79SLois Curfman McInnes 
3544a2ae208SSatish Balay #undef __FUNCT__
3554a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
356dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3578965ea79SLois Curfman McInnes {
358dfbe8321SBarry Smith   PetscErrorCode ierr;
35939ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3603a40ed3dSBarry Smith 
3613a40ed3dSBarry Smith   PetscFunctionBegin;
3623a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3633a40ed3dSBarry Smith   PetscFunctionReturn(0);
3648965ea79SLois Curfman McInnes }
3658965ea79SLois Curfman McInnes 
3668965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3678965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3688965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3693501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3708965ea79SLois Curfman McInnes    routine.
3718965ea79SLois Curfman McInnes */
3724a2ae208SSatish Balay #undef __FUNCT__
3734a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
374f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
3758965ea79SLois Curfman McInnes {
37639ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3776849ba73SBarry Smith   PetscErrorCode ierr;
378f4df32b1SMatthew Knepley   PetscInt       i,*owners = l->rowners;
37913f74950SBarry Smith   PetscInt       *nprocs,j,idx,nsends;
38013f74950SBarry Smith   PetscInt       nmax,*svalues,*starts,*owner,nrecvs;
38113f74950SBarry Smith   PetscInt       *rvalues,tag = A->tag,count,base,slen,*source;
38213f74950SBarry Smith   PetscInt       *lens,*lrows,*values;
38313f74950SBarry Smith   PetscMPIInt    n,imdex,rank = l->rank,size = l->size;
3848965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
3858965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
3868965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
38735d8aa7fSBarry Smith   PetscTruth     found;
3888965ea79SLois Curfman McInnes 
3893a40ed3dSBarry Smith   PetscFunctionBegin;
3908965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
39113f74950SBarry Smith   ierr  = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
39213f74950SBarry Smith   ierr  = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
39313f74950SBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/
3948965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3958965ea79SLois Curfman McInnes     idx = rows[i];
39635d8aa7fSBarry Smith     found = PETSC_FALSE;
3978965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3988965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
399c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
4008965ea79SLois Curfman McInnes       }
4018965ea79SLois Curfman McInnes     }
40229bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
4038965ea79SLois Curfman McInnes   }
404c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
4058965ea79SLois Curfman McInnes 
4068965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
407c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
4088965ea79SLois Curfman McInnes 
4098965ea79SLois Curfman McInnes   /* post receives:   */
41013f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
411b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
4128965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
41313f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
4148965ea79SLois Curfman McInnes   }
4158965ea79SLois Curfman McInnes 
4168965ea79SLois Curfman McInnes   /* do sends:
4178965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
4188965ea79SLois Curfman McInnes          the ith processor
4198965ea79SLois Curfman McInnes   */
42013f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
421b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
42213f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
4238965ea79SLois Curfman McInnes   starts[0]  = 0;
424c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
4258965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
4268965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
4278965ea79SLois Curfman McInnes   }
4288965ea79SLois Curfman McInnes 
4298965ea79SLois Curfman McInnes   starts[0] = 0;
430c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
4318965ea79SLois Curfman McInnes   count = 0;
4328965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
433c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
43413f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
4358965ea79SLois Curfman McInnes     }
4368965ea79SLois Curfman McInnes   }
437606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
4388965ea79SLois Curfman McInnes 
4398965ea79SLois Curfman McInnes   base = owners[rank];
4408965ea79SLois Curfman McInnes 
4418965ea79SLois Curfman McInnes   /*  wait on receives */
44213f74950SBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
4438965ea79SLois Curfman McInnes   source = lens + nrecvs;
4448965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
4458965ea79SLois Curfman McInnes   while (count) {
446ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4478965ea79SLois Curfman McInnes     /* unpack receives into our local space */
44813f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4498965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
4508965ea79SLois Curfman McInnes     lens[imdex]    = n;
4518965ea79SLois Curfman McInnes     slen += n;
4528965ea79SLois Curfman McInnes     count--;
4538965ea79SLois Curfman McInnes   }
454606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4558965ea79SLois Curfman McInnes 
4568965ea79SLois Curfman McInnes   /* move the data into the send scatter */
45713f74950SBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
4588965ea79SLois Curfman McInnes   count = 0;
4598965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4608965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4618965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4628965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4638965ea79SLois Curfman McInnes     }
4648965ea79SLois Curfman McInnes   }
465606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
466606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
467606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
468606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
4698965ea79SLois Curfman McInnes 
4708965ea79SLois Curfman McInnes   /* actually zap the local rows */
471f4df32b1SMatthew Knepley   ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr);
472606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4738965ea79SLois Curfman McInnes 
4748965ea79SLois Curfman McInnes   /* wait on sends */
4758965ea79SLois Curfman McInnes   if (nsends) {
476b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
477ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
478606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4798965ea79SLois Curfman McInnes   }
480606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
481606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4828965ea79SLois Curfman McInnes 
4833a40ed3dSBarry Smith   PetscFunctionReturn(0);
4848965ea79SLois Curfman McInnes }
4858965ea79SLois Curfman McInnes 
4864a2ae208SSatish Balay #undef __FUNCT__
4874a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
488dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4898965ea79SLois Curfman McInnes {
49039ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
491dfbe8321SBarry Smith   PetscErrorCode ierr;
492c456f294SBarry Smith 
4933a40ed3dSBarry Smith   PetscFunctionBegin;
49443a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
49543a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
49644cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4973a40ed3dSBarry Smith   PetscFunctionReturn(0);
4988965ea79SLois Curfman McInnes }
4998965ea79SLois Curfman McInnes 
5004a2ae208SSatish Balay #undef __FUNCT__
5014a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
502dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
5038965ea79SLois Curfman McInnes {
50439ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
505dfbe8321SBarry Smith   PetscErrorCode ierr;
506c456f294SBarry Smith 
5073a40ed3dSBarry Smith   PetscFunctionBegin;
50843a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
50943a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
51044cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
5113a40ed3dSBarry Smith   PetscFunctionReturn(0);
5128965ea79SLois Curfman McInnes }
5138965ea79SLois Curfman McInnes 
5144a2ae208SSatish Balay #undef __FUNCT__
5154a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
516dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
517096963f5SLois Curfman McInnes {
518096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
519dfbe8321SBarry Smith   PetscErrorCode ierr;
52087828ca2SBarry Smith   PetscScalar    zero = 0.0;
521096963f5SLois Curfman McInnes 
5223a40ed3dSBarry Smith   PetscFunctionBegin;
5232dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
5247c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
525537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
526537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
5273a40ed3dSBarry Smith   PetscFunctionReturn(0);
528096963f5SLois Curfman McInnes }
529096963f5SLois Curfman McInnes 
5304a2ae208SSatish Balay #undef __FUNCT__
5314a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
532dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
533096963f5SLois Curfman McInnes {
534096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
535dfbe8321SBarry Smith   PetscErrorCode ierr;
536096963f5SLois Curfman McInnes 
5373a40ed3dSBarry Smith   PetscFunctionBegin;
5383501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
5397c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
540537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
541537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
5423a40ed3dSBarry Smith   PetscFunctionReturn(0);
543096963f5SLois Curfman McInnes }
544096963f5SLois Curfman McInnes 
5454a2ae208SSatish Balay #undef __FUNCT__
5464a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
547dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5488965ea79SLois Curfman McInnes {
54939ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
550096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
551dfbe8321SBarry Smith   PetscErrorCode ierr;
55213f74950SBarry Smith   PetscInt       len,i,n,m = A->m,radd;
55387828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
554ed3cc1f0SBarry Smith 
5553a40ed3dSBarry Smith   PetscFunctionBegin;
5562dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5571ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
558096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
559273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
560273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
5617ddc982cSLois Curfman McInnes   radd = a->rstart*m;
56244cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
563096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
564096963f5SLois Curfman McInnes   }
5651ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5663a40ed3dSBarry Smith   PetscFunctionReturn(0);
5678965ea79SLois Curfman McInnes }
5688965ea79SLois Curfman McInnes 
5694a2ae208SSatish Balay #undef __FUNCT__
5704a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
571dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5728965ea79SLois Curfman McInnes {
5733501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
574dfbe8321SBarry Smith   PetscErrorCode ierr;
575ed3cc1f0SBarry Smith 
5763a40ed3dSBarry Smith   PetscFunctionBegin;
57794d884c6SBarry Smith 
578aa482453SBarry Smith #if defined(PETSC_USE_LOG)
57977431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N);
5808965ea79SLois Curfman McInnes #endif
5818798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
582606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
5833501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
5843501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
5853501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
586622d7880SLois Curfman McInnes   if (mdn->factor) {
587606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
588606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
589606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
590606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
591622d7880SLois Curfman McInnes   }
592606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
593901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
594901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5953a40ed3dSBarry Smith   PetscFunctionReturn(0);
5968965ea79SLois Curfman McInnes }
59739ddd567SLois Curfman McInnes 
5984a2ae208SSatish Balay #undef __FUNCT__
5994a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
6006849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
6018965ea79SLois Curfman McInnes {
60239ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
603dfbe8321SBarry Smith   PetscErrorCode ierr;
6047056b6fcSBarry Smith 
6053a40ed3dSBarry Smith   PetscFunctionBegin;
60639ddd567SLois Curfman McInnes   if (mdn->size == 1) {
60739ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6088965ea79SLois Curfman McInnes   }
60929bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
6103a40ed3dSBarry Smith   PetscFunctionReturn(0);
6118965ea79SLois Curfman McInnes }
6128965ea79SLois Curfman McInnes 
6134a2ae208SSatish Balay #undef __FUNCT__
6144a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
6156849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6168965ea79SLois Curfman McInnes {
61739ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
618dfbe8321SBarry Smith   PetscErrorCode    ierr;
61932dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
620b0a32e0cSBarry Smith   PetscViewerType   vtype;
62132077d6dSBarry Smith   PetscTruth        iascii,isdraw;
622b0a32e0cSBarry Smith   PetscViewer       sviewer;
623f3ef73ceSBarry Smith   PetscViewerFormat format;
6248965ea79SLois Curfman McInnes 
6253a40ed3dSBarry Smith   PetscFunctionBegin;
62632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
627fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
62832077d6dSBarry Smith   if (iascii) {
629b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
630b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
631456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6324e220ebcSLois Curfman McInnes       MatInfo info;
633888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
63477431f27SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m,
63577431f27SBarry Smith                    (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr);
636b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6373501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
6383a40ed3dSBarry Smith       PetscFunctionReturn(0);
639fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6403a40ed3dSBarry Smith       PetscFunctionReturn(0);
6418965ea79SLois Curfman McInnes     }
642f1af5d2fSBarry Smith   } else if (isdraw) {
643b0a32e0cSBarry Smith     PetscDraw       draw;
644f1af5d2fSBarry Smith     PetscTruth isnull;
645f1af5d2fSBarry Smith 
646b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
647b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
648f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
649f1af5d2fSBarry Smith   }
65077ed5343SBarry Smith 
6518965ea79SLois Curfman McInnes   if (size == 1) {
65239ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6533a40ed3dSBarry Smith   } else {
6548965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
6558965ea79SLois Curfman McInnes     Mat         A;
65613f74950SBarry Smith     PetscInt    M = mat->M,N = mat->N,m,row,i,nz;
657ba8c8a56SBarry Smith     PetscInt    *cols;
658ba8c8a56SBarry Smith     PetscScalar *vals;
6598965ea79SLois Curfman McInnes 
660f69a0ea3SMatthew Knepley     ierr = MatCreate(mat->comm,&A);CHKERRQ(ierr);
6618965ea79SLois Curfman McInnes     if (!rank) {
662f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
6633a40ed3dSBarry Smith     } else {
664f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
6658965ea79SLois Curfman McInnes     }
666be5d1d56SKris Buschelman     /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */
667878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
668878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
66952e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
6708965ea79SLois Curfman McInnes 
67139ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
67239ddd567SLois Curfman McInnes        but it's quick for now */
67351022da4SBarry Smith     A->insertmode = INSERT_VALUES;
674273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
6758965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
676ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
677ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
678ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
67939ddd567SLois Curfman McInnes       row++;
6808965ea79SLois Curfman McInnes     }
6818965ea79SLois Curfman McInnes 
6826d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6836d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
684b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
685b9b97703SBarry Smith     if (!rank) {
6866831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
6878965ea79SLois Curfman McInnes     }
688b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
689b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6908965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
6918965ea79SLois Curfman McInnes   }
6923a40ed3dSBarry Smith   PetscFunctionReturn(0);
6938965ea79SLois Curfman McInnes }
6948965ea79SLois Curfman McInnes 
6954a2ae208SSatish Balay #undef __FUNCT__
6964a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
697dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
6988965ea79SLois Curfman McInnes {
699dfbe8321SBarry Smith   PetscErrorCode ierr;
70032077d6dSBarry Smith   PetscTruth     iascii,isbinary,isdraw,issocket;
7018965ea79SLois Curfman McInnes 
702433994e6SBarry Smith   PetscFunctionBegin;
7030f5bd95cSBarry Smith 
70432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
705fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
706b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
707fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
7080f5bd95cSBarry Smith 
70932077d6dSBarry Smith   if (iascii || issocket || isdraw) {
710f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7110f5bd95cSBarry Smith   } else if (isbinary) {
7123a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
7135cd90555SBarry Smith   } else {
714958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
7158965ea79SLois Curfman McInnes   }
7163a40ed3dSBarry Smith   PetscFunctionReturn(0);
7178965ea79SLois Curfman McInnes }
7188965ea79SLois Curfman McInnes 
7194a2ae208SSatish Balay #undef __FUNCT__
7204a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
721dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7228965ea79SLois Curfman McInnes {
7233501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7243501a2bdSLois Curfman McInnes   Mat            mdn = mat->A;
725dfbe8321SBarry Smith   PetscErrorCode ierr;
726329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7278965ea79SLois Curfman McInnes 
7283a40ed3dSBarry Smith   PetscFunctionBegin;
729273d9f13SBarry Smith   info->rows_global    = (double)A->M;
730273d9f13SBarry Smith   info->columns_global = (double)A->N;
731273d9f13SBarry Smith   info->rows_local     = (double)A->m;
732273d9f13SBarry Smith   info->columns_local  = (double)A->N;
7334e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
7344e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7354e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7364e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7378965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7384e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7394e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7404e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7414e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7424e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7438965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
744d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
7454e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7464e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7474e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7484e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7494e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7508965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
751d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
7524e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7534e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7544e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7554e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7564e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7578965ea79SLois Curfman McInnes   }
7584e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
7594e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
7604e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
7613a40ed3dSBarry Smith   PetscFunctionReturn(0);
7628965ea79SLois Curfman McInnes }
7638965ea79SLois Curfman McInnes 
7644a2ae208SSatish Balay #undef __FUNCT__
7654a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
766dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op)
7678965ea79SLois Curfman McInnes {
76839ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
769dfbe8321SBarry Smith   PetscErrorCode ierr;
7708965ea79SLois Curfman McInnes 
7713a40ed3dSBarry Smith   PetscFunctionBegin;
77212c028f9SKris Buschelman   switch (op) {
77312c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
77412c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
77512c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
77612c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
77712c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
77812c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
779273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
78012c028f9SKris Buschelman     break;
78112c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
782273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
783273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
78412c028f9SKris Buschelman     break;
78512c028f9SKris Buschelman   case MAT_ROWS_SORTED:
78612c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
78712c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
78812c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
78963ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatSetOption_MPIDense:Option ignored\n"));CHKERRQ(ierr);
79012c028f9SKris Buschelman     break;
79112c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
792273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
793273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
79412c028f9SKris Buschelman     break;
79512c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
796273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
79712c028f9SKris Buschelman     break;
79812c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
79929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
80077e54ba9SKris Buschelman   case MAT_SYMMETRIC:
80177e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8029a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
8039a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
8049a4540c5SBarry Smith   case MAT_HERMITIAN:
8059a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
8069a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
8079a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
80877e54ba9SKris Buschelman     break;
80912c028f9SKris Buschelman   default:
81029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
8113a40ed3dSBarry Smith   }
8123a40ed3dSBarry Smith   PetscFunctionReturn(0);
8138965ea79SLois Curfman McInnes }
8148965ea79SLois Curfman McInnes 
8158965ea79SLois Curfman McInnes 
8164a2ae208SSatish Balay #undef __FUNCT__
8174a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
818dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8195b2fa520SLois Curfman McInnes {
8205b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8215b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
82287828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
823dfbe8321SBarry Smith   PetscErrorCode ierr;
82413f74950SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
8255b2fa520SLois Curfman McInnes 
8265b2fa520SLois Curfman McInnes   PetscFunctionBegin;
82772d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8285b2fa520SLois Curfman McInnes   if (ll) {
82972d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
830*175be7b4SMatthew Knepley     if (s2a != s2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8311ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8325b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8335b2fa520SLois Curfman McInnes       x = l[i];
8345b2fa520SLois Curfman McInnes       v = mat->v + i;
8355b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8365b2fa520SLois Curfman McInnes     }
8371ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
838efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8395b2fa520SLois Curfman McInnes   }
8405b2fa520SLois Curfman McInnes   if (rr) {
841*175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
842*175be7b4SMatthew Knepley     if (s3a != s3) SETERRQ2(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
8435b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8445b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8451ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8465b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8475b2fa520SLois Curfman McInnes       x = r[i];
8485b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8495b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
8505b2fa520SLois Curfman McInnes     }
8511ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
852efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8535b2fa520SLois Curfman McInnes   }
8545b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8555b2fa520SLois Curfman McInnes }
8565b2fa520SLois Curfman McInnes 
8574a2ae208SSatish Balay #undef __FUNCT__
8584a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
859dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
860096963f5SLois Curfman McInnes {
8613501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8623501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
863dfbe8321SBarry Smith   PetscErrorCode ierr;
86413f74950SBarry Smith   PetscInt       i,j;
865329f5518SBarry Smith   PetscReal      sum = 0.0;
86687828ca2SBarry Smith   PetscScalar    *v = mat->v;
8673501a2bdSLois Curfman McInnes 
8683a40ed3dSBarry Smith   PetscFunctionBegin;
8693501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
870064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
8713501a2bdSLois Curfman McInnes   } else {
8723501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
873273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
874aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
875329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
8763501a2bdSLois Curfman McInnes #else
8773501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
8783501a2bdSLois Curfman McInnes #endif
8793501a2bdSLois Curfman McInnes       }
880064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
881064f8208SBarry Smith       *nrm = sqrt(*nrm);
882efee365bSSatish Balay       ierr = PetscLogFlops(2*mdn->A->n*mdn->A->m);CHKERRQ(ierr);
8833a40ed3dSBarry Smith     } else if (type == NORM_1) {
884329f5518SBarry Smith       PetscReal *tmp,*tmp2;
885b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
886273d9f13SBarry Smith       tmp2 = tmp + A->N;
887273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
888064f8208SBarry Smith       *nrm = 0.0;
8893501a2bdSLois Curfman McInnes       v = mat->v;
890273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
891273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
89267e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8933501a2bdSLois Curfman McInnes         }
8943501a2bdSLois Curfman McInnes       }
895d7d1e502SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
896273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
897064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8983501a2bdSLois Curfman McInnes       }
899606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
900efee365bSSatish Balay       ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
9013a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
902329f5518SBarry Smith       PetscReal ntemp;
9033501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
904064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
9053a40ed3dSBarry Smith     } else {
90629bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
9073501a2bdSLois Curfman McInnes     }
9083501a2bdSLois Curfman McInnes   }
9093a40ed3dSBarry Smith   PetscFunctionReturn(0);
9103501a2bdSLois Curfman McInnes }
9113501a2bdSLois Curfman McInnes 
9124a2ae208SSatish Balay #undef __FUNCT__
9134a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
914dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout)
9153501a2bdSLois Curfman McInnes {
9163501a2bdSLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
9173501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
9183501a2bdSLois Curfman McInnes   Mat            B;
91913f74950SBarry Smith   PetscInt       M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
9206849ba73SBarry Smith   PetscErrorCode ierr;
92113f74950SBarry Smith   PetscInt       j,i;
92287828ca2SBarry Smith   PetscScalar    *v;
9233501a2bdSLois Curfman McInnes 
9243a40ed3dSBarry Smith   PetscFunctionBegin;
9257c922b88SBarry Smith   if (!matout && M != N) {
92629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
9277056b6fcSBarry Smith   }
928f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
929f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr);
930878740d9SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
931878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
9323501a2bdSLois Curfman McInnes 
933273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
93413f74950SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
9353501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
9363501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
9373501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9383501a2bdSLois Curfman McInnes     v   += m;
9393501a2bdSLois Curfman McInnes   }
940606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9416d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9426d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9437c922b88SBarry Smith   if (matout) {
9443501a2bdSLois Curfman McInnes     *matout = B;
9453501a2bdSLois Curfman McInnes   } else {
946273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
9473501a2bdSLois Curfman McInnes   }
9483a40ed3dSBarry Smith   PetscFunctionReturn(0);
949096963f5SLois Curfman McInnes }
950096963f5SLois Curfman McInnes 
951d9eff348SSatish Balay #include "petscblaslapack.h"
9524a2ae208SSatish Balay #undef __FUNCT__
9534a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
954f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha)
95544cd7ae7SLois Curfman McInnes {
95644cd7ae7SLois Curfman McInnes   Mat_MPIDense   *A = (Mat_MPIDense*)inA->data;
95744cd7ae7SLois Curfman McInnes   Mat_SeqDense   *a = (Mat_SeqDense*)A->A->data;
958f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
9594ce68768SBarry Smith   PetscBLASInt   one = 1,nz = (PetscBLASInt)inA->m*inA->N;
960efee365bSSatish Balay   PetscErrorCode ierr;
96144cd7ae7SLois Curfman McInnes 
9623a40ed3dSBarry Smith   PetscFunctionBegin;
963f4df32b1SMatthew Knepley   BLASscal_(&nz,&oalpha,a->v,&one);
964efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
9653a40ed3dSBarry Smith   PetscFunctionReturn(0);
96644cd7ae7SLois Curfman McInnes }
96744cd7ae7SLois Curfman McInnes 
9686849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
9698965ea79SLois Curfman McInnes 
9704a2ae208SSatish Balay #undef __FUNCT__
9714a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
972dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
973273d9f13SBarry Smith {
974dfbe8321SBarry Smith   PetscErrorCode ierr;
975273d9f13SBarry Smith 
976273d9f13SBarry Smith   PetscFunctionBegin;
977273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
978273d9f13SBarry Smith   PetscFunctionReturn(0);
979273d9f13SBarry Smith }
980273d9f13SBarry Smith 
9818965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
98209dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
98309dc0095SBarry Smith        MatGetRow_MPIDense,
98409dc0095SBarry Smith        MatRestoreRow_MPIDense,
98509dc0095SBarry Smith        MatMult_MPIDense,
98697304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
9877c922b88SBarry Smith        MatMultTranspose_MPIDense,
9887c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9898965ea79SLois Curfman McInnes        0,
99009dc0095SBarry Smith        0,
99109dc0095SBarry Smith        0,
99297304618SKris Buschelman /*10*/ 0,
99309dc0095SBarry Smith        0,
99409dc0095SBarry Smith        0,
99509dc0095SBarry Smith        0,
99609dc0095SBarry Smith        MatTranspose_MPIDense,
99797304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
99897304618SKris Buschelman        0,
99909dc0095SBarry Smith        MatGetDiagonal_MPIDense,
10005b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
100109dc0095SBarry Smith        MatNorm_MPIDense,
100297304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
100309dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
100409dc0095SBarry Smith        0,
100509dc0095SBarry Smith        MatSetOption_MPIDense,
100609dc0095SBarry Smith        MatZeroEntries_MPIDense,
100797304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
10080ad19415SHong Zhang        MatLUFactorSymbolic_MPIDense,
10090ad19415SHong Zhang        MatLUFactorNumeric_MPIDense,
10100ad19415SHong Zhang        MatCholeskyFactorSymbolic_MPIDense,
10110ad19415SHong Zhang        MatCholeskyFactorNumeric_MPIDense,
101297304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
1013273d9f13SBarry Smith        0,
101409dc0095SBarry Smith        0,
101509dc0095SBarry Smith        MatGetArray_MPIDense,
101609dc0095SBarry Smith        MatRestoreArray_MPIDense,
101797304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
101809dc0095SBarry Smith        0,
101909dc0095SBarry Smith        0,
102009dc0095SBarry Smith        0,
102109dc0095SBarry Smith        0,
102297304618SKris Buschelman /*40*/ 0,
10232ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
102409dc0095SBarry Smith        0,
102509dc0095SBarry Smith        MatGetValues_MPIDense,
102609dc0095SBarry Smith        0,
102797304618SKris Buschelman /*45*/ 0,
102809dc0095SBarry Smith        MatScale_MPIDense,
102909dc0095SBarry Smith        0,
103009dc0095SBarry Smith        0,
103109dc0095SBarry Smith        0,
1032521d7252SBarry Smith /*50*/ 0,
103309dc0095SBarry Smith        0,
103409dc0095SBarry Smith        0,
103509dc0095SBarry Smith        0,
103609dc0095SBarry Smith        0,
103797304618SKris Buschelman /*55*/ 0,
103809dc0095SBarry Smith        0,
103909dc0095SBarry Smith        0,
104009dc0095SBarry Smith        0,
104109dc0095SBarry Smith        0,
104297304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
1043b9b97703SBarry Smith        MatDestroy_MPIDense,
1044b9b97703SBarry Smith        MatView_MPIDense,
104597304618SKris Buschelman        MatGetPetscMaps_Petsc,
104697304618SKris Buschelman        0,
104797304618SKris Buschelman /*65*/ 0,
104897304618SKris Buschelman        0,
104997304618SKris Buschelman        0,
105097304618SKris Buschelman        0,
105197304618SKris Buschelman        0,
105297304618SKris Buschelman /*70*/ 0,
105397304618SKris Buschelman        0,
105497304618SKris Buschelman        0,
105597304618SKris Buschelman        0,
105697304618SKris Buschelman        0,
105797304618SKris Buschelman /*75*/ 0,
105897304618SKris Buschelman        0,
105997304618SKris Buschelman        0,
106097304618SKris Buschelman        0,
106197304618SKris Buschelman        0,
106297304618SKris Buschelman /*80*/ 0,
106397304618SKris Buschelman        0,
106497304618SKris Buschelman        0,
106597304618SKris Buschelman        0,
1066865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense,
1067865e5f61SKris Buschelman        0,
1068865e5f61SKris Buschelman        0,
1069865e5f61SKris Buschelman        0,
1070865e5f61SKris Buschelman        0,
1071865e5f61SKris Buschelman        0,
1072865e5f61SKris Buschelman /*90*/ 0,
1073865e5f61SKris Buschelman        0,
1074865e5f61SKris Buschelman        0,
1075865e5f61SKris Buschelman        0,
1076865e5f61SKris Buschelman        0,
1077865e5f61SKris Buschelman /*95*/ 0,
1078865e5f61SKris Buschelman        0,
1079865e5f61SKris Buschelman        0,
1080865e5f61SKris Buschelman        0};
10818965ea79SLois Curfman McInnes 
1082273d9f13SBarry Smith EXTERN_C_BEGIN
10834a2ae208SSatish Balay #undef __FUNCT__
1084a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
1085be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1086a23d5eceSKris Buschelman {
1087a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1088dfbe8321SBarry Smith   PetscErrorCode ierr;
1089a23d5eceSKris Buschelman 
1090a23d5eceSKris Buschelman   PetscFunctionBegin;
1091a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1092a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1093a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1094a23d5eceSKris Buschelman 
1095a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
1096f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1097f69a0ea3SMatthew Knepley   ierr = MatSetSizes(a->A,mat->m,mat->N,mat->m,mat->N);CHKERRQ(ierr);
10985c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
10995c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
110052e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1101a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1102a23d5eceSKris Buschelman }
1103a23d5eceSKris Buschelman EXTERN_C_END
1104a23d5eceSKris Buschelman 
11050bad9183SKris Buschelman /*MC
1106fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
11070bad9183SKris Buschelman 
11080bad9183SKris Buschelman    Options Database Keys:
11090bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
11100bad9183SKris Buschelman 
11110bad9183SKris Buschelman   Level: beginner
11120bad9183SKris Buschelman 
11130bad9183SKris Buschelman .seealso: MatCreateMPIDense
11140bad9183SKris Buschelman M*/
11150bad9183SKris Buschelman 
1116a23d5eceSKris Buschelman EXTERN_C_BEGIN
1117a23d5eceSKris Buschelman #undef __FUNCT__
11184a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1119be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat)
1120273d9f13SBarry Smith {
1121273d9f13SBarry Smith   Mat_MPIDense   *a;
1122dfbe8321SBarry Smith   PetscErrorCode ierr;
112313f74950SBarry Smith   PetscInt       i;
1124273d9f13SBarry Smith 
1125273d9f13SBarry Smith   PetscFunctionBegin;
1126b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1127b0a32e0cSBarry Smith   mat->data         = (void*)a;
1128273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1129273d9f13SBarry Smith   mat->factor       = 0;
1130273d9f13SBarry Smith   mat->mapping      = 0;
1131273d9f13SBarry Smith 
1132273d9f13SBarry Smith   a->factor       = 0;
1133273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1134273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
1135273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
1136273d9f13SBarry Smith 
1137273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
1138273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
1139273d9f13SBarry Smith   a->nvec = mat->n;
1140273d9f13SBarry Smith 
1141273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1142273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
11438a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
11448a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
1145273d9f13SBarry Smith 
1146273d9f13SBarry Smith   /* build local table of row and column ownerships */
114713f74950SBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr);
1148273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
114952e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr);
115013f74950SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr);
1151273d9f13SBarry Smith   a->rowners[0] = 0;
1152273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1153273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
1154273d9f13SBarry Smith   }
1155273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
1156273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
115713f74950SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr);
1158273d9f13SBarry Smith   a->cowners[0] = 0;
1159273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1160273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1161273d9f13SBarry Smith   }
1162273d9f13SBarry Smith 
1163273d9f13SBarry Smith   /* build cache for off array entries formed */
1164273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1165273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1166273d9f13SBarry Smith 
1167273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1168273d9f13SBarry Smith   a->lvec        = 0;
1169273d9f13SBarry Smith   a->Mvctx       = 0;
1170273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1171273d9f13SBarry Smith 
1172273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1173273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1174273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1175a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1176a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1177a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1178273d9f13SBarry Smith   PetscFunctionReturn(0);
1179273d9f13SBarry Smith }
1180273d9f13SBarry Smith EXTERN_C_END
1181273d9f13SBarry Smith 
1182209238afSKris Buschelman /*MC
1183002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1184209238afSKris Buschelman 
1185209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1186209238afSKris Buschelman    and MATMPIDENSE otherwise.
1187209238afSKris Buschelman 
1188209238afSKris Buschelman    Options Database Keys:
1189209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1190209238afSKris Buschelman 
1191209238afSKris Buschelman   Level: beginner
1192209238afSKris Buschelman 
1193209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1194209238afSKris Buschelman M*/
1195209238afSKris Buschelman 
1196209238afSKris Buschelman EXTERN_C_BEGIN
1197209238afSKris Buschelman #undef __FUNCT__
1198209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1199be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A)
1200dfbe8321SBarry Smith {
12016849ba73SBarry Smith   PetscErrorCode ierr;
120213f74950SBarry Smith   PetscMPIInt    size;
1203209238afSKris Buschelman 
1204209238afSKris Buschelman   PetscFunctionBegin;
1205209238afSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr);
1206209238afSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1207209238afSKris Buschelman   if (size == 1) {
1208209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1209209238afSKris Buschelman   } else {
1210209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1211209238afSKris Buschelman   }
1212209238afSKris Buschelman   PetscFunctionReturn(0);
1213209238afSKris Buschelman }
1214209238afSKris Buschelman EXTERN_C_END
1215209238afSKris Buschelman 
12164a2ae208SSatish Balay #undef __FUNCT__
12174a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1218273d9f13SBarry Smith /*@C
1219273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1220273d9f13SBarry Smith 
1221273d9f13SBarry Smith    Not collective
1222273d9f13SBarry Smith 
1223273d9f13SBarry Smith    Input Parameters:
1224273d9f13SBarry Smith .  A - the matrix
1225273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1226273d9f13SBarry Smith    to control all matrix memory allocation.
1227273d9f13SBarry Smith 
1228273d9f13SBarry Smith    Notes:
1229273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1230273d9f13SBarry Smith    storage by columns.
1231273d9f13SBarry Smith 
1232273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1233273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1234273d9f13SBarry Smith    set data=PETSC_NULL.
1235273d9f13SBarry Smith 
1236273d9f13SBarry Smith    Level: intermediate
1237273d9f13SBarry Smith 
1238273d9f13SBarry Smith .keywords: matrix,dense, parallel
1239273d9f13SBarry Smith 
1240273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1241273d9f13SBarry Smith @*/
1242be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1243273d9f13SBarry Smith {
1244dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar *);
1245273d9f13SBarry Smith 
1246273d9f13SBarry Smith   PetscFunctionBegin;
1247565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1248a23d5eceSKris Buschelman   if (f) {
1249a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1250a23d5eceSKris Buschelman   }
1251273d9f13SBarry Smith   PetscFunctionReturn(0);
1252273d9f13SBarry Smith }
1253273d9f13SBarry Smith 
12544a2ae208SSatish Balay #undef __FUNCT__
12554a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
12568965ea79SLois Curfman McInnes /*@C
125739ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
12588965ea79SLois Curfman McInnes 
1259db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1260db81eaa0SLois Curfman McInnes 
12618965ea79SLois Curfman McInnes    Input Parameters:
1262db81eaa0SLois Curfman McInnes +  comm - MPI communicator
12638965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1264db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
12658965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1266db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
12677f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1268dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
12698965ea79SLois Curfman McInnes 
12708965ea79SLois Curfman McInnes    Output Parameter:
1271477f1c0bSLois Curfman McInnes .  A - the matrix
12728965ea79SLois Curfman McInnes 
1273b259b22eSLois Curfman McInnes    Notes:
127439ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
127539ddd567SLois Curfman McInnes    storage by columns.
12768965ea79SLois Curfman McInnes 
127718f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
127818f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
12797f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
128018f449edSLois Curfman McInnes 
12818965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
12828965ea79SLois Curfman McInnes    (possibly both).
12838965ea79SLois Curfman McInnes 
1284027ccd11SLois Curfman McInnes    Level: intermediate
1285027ccd11SLois Curfman McInnes 
128639ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
12878965ea79SLois Curfman McInnes 
128839ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
12898965ea79SLois Curfman McInnes @*/
1290be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
12918965ea79SLois Curfman McInnes {
12926849ba73SBarry Smith   PetscErrorCode ierr;
129313f74950SBarry Smith   PetscMPIInt    size;
12948965ea79SLois Curfman McInnes 
12953a40ed3dSBarry Smith   PetscFunctionBegin;
1296f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1297f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1298273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1299273d9f13SBarry Smith   if (size > 1) {
1300273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1301273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1302273d9f13SBarry Smith   } else {
1303273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1304273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
13058c469469SLois Curfman McInnes   }
13063a40ed3dSBarry Smith   PetscFunctionReturn(0);
13078965ea79SLois Curfman McInnes }
13088965ea79SLois Curfman McInnes 
13094a2ae208SSatish Balay #undef __FUNCT__
13104a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
13116849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
13128965ea79SLois Curfman McInnes {
13138965ea79SLois Curfman McInnes   Mat            mat;
13143501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1315dfbe8321SBarry Smith   PetscErrorCode ierr;
13168965ea79SLois Curfman McInnes 
13173a40ed3dSBarry Smith   PetscFunctionBegin;
13188965ea79SLois Curfman McInnes   *newmat       = 0;
1319f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&mat);CHKERRQ(ierr);
1320f69a0ea3SMatthew Knepley   ierr = MatSetSizes(mat,A->m,A->n,A->M,A->N);CHKERRQ(ierr);
1321be5d1d56SKris Buschelman   ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr);
1322834f8fabSBarry Smith   a                 = (Mat_MPIDense*)mat->data;
1323e04c1aa4SHong Zhang   ierr              = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
13243501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1325c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1326273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
13278965ea79SLois Curfman McInnes 
13288965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
13298965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
13308965ea79SLois Curfman McInnes   a->size         = oldmat->size;
13318965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1332e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1333b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
13343782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1335e04c1aa4SHong Zhang 
133652e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr);
133713f74950SBarry Smith   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
13388798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
13398965ea79SLois Curfman McInnes 
1340329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
13415609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
134252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
13438965ea79SLois Curfman McInnes   *newmat = mat;
13443a40ed3dSBarry Smith   PetscFunctionReturn(0);
13458965ea79SLois Curfman McInnes }
13468965ea79SLois Curfman McInnes 
1347e090d566SSatish Balay #include "petscsys.h"
13488965ea79SLois Curfman McInnes 
13494a2ae208SSatish Balay #undef __FUNCT__
13504a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
1351f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat)
135290ace30eSBarry Smith {
13536849ba73SBarry Smith   PetscErrorCode ierr;
135432dcc486SBarry Smith   PetscMPIInt    rank,size;
135513f74950SBarry Smith   PetscInt       *rowners,i,m,nz,j;
135687828ca2SBarry Smith   PetscScalar    *array,*vals,*vals_ptr;
135790ace30eSBarry Smith   MPI_Status     status;
135890ace30eSBarry Smith 
13593a40ed3dSBarry Smith   PetscFunctionBegin;
1360d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1361d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
136290ace30eSBarry Smith 
136390ace30eSBarry Smith   /* determine ownership of all rows */
136490ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
136513f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
136613f74950SBarry Smith   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
136790ace30eSBarry Smith   rowners[0] = 0;
136890ace30eSBarry Smith   for (i=2; i<=size; i++) {
136990ace30eSBarry Smith     rowners[i] += rowners[i-1];
137090ace30eSBarry Smith   }
137190ace30eSBarry Smith 
1372f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1373f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1374be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1375878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
137690ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
137790ace30eSBarry Smith 
137890ace30eSBarry Smith   if (!rank) {
137987828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
138090ace30eSBarry Smith 
138190ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
13820752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
138390ace30eSBarry Smith 
138490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
138590ace30eSBarry Smith     vals_ptr = vals;
138690ace30eSBarry Smith     for (i=0; i<m; i++) {
138790ace30eSBarry Smith       for (j=0; j<N; j++) {
138890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
138990ace30eSBarry Smith       }
139090ace30eSBarry Smith     }
139190ace30eSBarry Smith 
139290ace30eSBarry Smith     /* read in other processors and ship out */
139390ace30eSBarry Smith     for (i=1; i<size; i++) {
139490ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
13950752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1396ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
139790ace30eSBarry Smith     }
13983a40ed3dSBarry Smith   } else {
139990ace30eSBarry Smith     /* receive numeric values */
140087828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
140190ace30eSBarry Smith 
140290ace30eSBarry Smith     /* receive message of values*/
1403ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
140490ace30eSBarry Smith 
140590ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
140690ace30eSBarry Smith     vals_ptr = vals;
140790ace30eSBarry Smith     for (i=0; i<m; i++) {
140890ace30eSBarry Smith       for (j=0; j<N; j++) {
140990ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
141090ace30eSBarry Smith       }
141190ace30eSBarry Smith     }
141290ace30eSBarry Smith   }
1413606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1414606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
14156d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14166d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14173a40ed3dSBarry Smith   PetscFunctionReturn(0);
141890ace30eSBarry Smith }
141990ace30eSBarry Smith 
14204a2ae208SSatish Balay #undef __FUNCT__
14214a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
1422f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat)
14238965ea79SLois Curfman McInnes {
14248965ea79SLois Curfman McInnes   Mat            A;
142587828ca2SBarry Smith   PetscScalar    *vals,*svals;
142619bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
14278965ea79SLois Curfman McInnes   MPI_Status     status;
142813f74950SBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
142913f74950SBarry Smith   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
143013f74950SBarry Smith   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
143113f74950SBarry Smith   PetscInt       i,nz,j,rstart,rend;
143213f74950SBarry Smith   int            fd;
14336849ba73SBarry Smith   PetscErrorCode ierr;
14348965ea79SLois Curfman McInnes 
14353a40ed3dSBarry Smith   PetscFunctionBegin;
1436d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1437d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
14388965ea79SLois Curfman McInnes   if (!rank) {
1439b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
14400752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1441552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
14428965ea79SLois Curfman McInnes   }
14438965ea79SLois Curfman McInnes 
144413f74950SBarry Smith   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
144590ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
144690ace30eSBarry Smith 
144790ace30eSBarry Smith   /*
144890ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
144990ace30eSBarry Smith   */
145090ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1451be5d1d56SKris Buschelman     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr);
14523a40ed3dSBarry Smith     PetscFunctionReturn(0);
145390ace30eSBarry Smith   }
145490ace30eSBarry Smith 
14558965ea79SLois Curfman McInnes   /* determine ownership of all rows */
14568965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
145713f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
1458ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
14598965ea79SLois Curfman McInnes   rowners[0] = 0;
14608965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
14618965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
14628965ea79SLois Curfman McInnes   }
14638965ea79SLois Curfman McInnes   rstart = rowners[rank];
14648965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
14658965ea79SLois Curfman McInnes 
14668965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
146713f74950SBarry Smith   ierr    = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr);
14688965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
14698965ea79SLois Curfman McInnes   if (!rank) {
147013f74950SBarry Smith     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
14710752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
147213f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
14738965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
14741466f1e1SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
1475606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1476ca161407SBarry Smith   } else {
14771466f1e1SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
14788965ea79SLois Curfman McInnes   }
14798965ea79SLois Curfman McInnes 
14808965ea79SLois Curfman McInnes   if (!rank) {
14818965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
148213f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
148313f74950SBarry Smith     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
14848965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14858965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
14868965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
14878965ea79SLois Curfman McInnes       }
14888965ea79SLois Curfman McInnes     }
1489606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
14908965ea79SLois Curfman McInnes 
14918965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
14928965ea79SLois Curfman McInnes     maxnz = 0;
14938965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14940452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
14958965ea79SLois Curfman McInnes     }
149613f74950SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
14978965ea79SLois Curfman McInnes 
14988965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
14998965ea79SLois Curfman McInnes     nz = procsnz[0];
150013f74950SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
15010752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
15028965ea79SLois Curfman McInnes 
15038965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
15048965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
15058965ea79SLois Curfman McInnes       nz   = procsnz[i];
15060752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
150713f74950SBarry Smith       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
15088965ea79SLois Curfman McInnes     }
1509606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
15103a40ed3dSBarry Smith   } else {
15118965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
15128965ea79SLois Curfman McInnes     nz = 0;
15138965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15148965ea79SLois Curfman McInnes       nz += ourlens[i];
15158965ea79SLois Curfman McInnes     }
151613f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
15178965ea79SLois Curfman McInnes 
15188965ea79SLois Curfman McInnes     /* receive message of column indices*/
151913f74950SBarry Smith     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
152013f74950SBarry Smith     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
152129bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
15228965ea79SLois Curfman McInnes   }
15238965ea79SLois Curfman McInnes 
15248965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
152513f74950SBarry Smith   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
15268965ea79SLois Curfman McInnes   jj = 0;
15278965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
15288965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
15298965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
15308965ea79SLois Curfman McInnes       jj++;
15318965ea79SLois Curfman McInnes     }
15328965ea79SLois Curfman McInnes   }
15338965ea79SLois Curfman McInnes 
15348965ea79SLois Curfman McInnes   /* create our matrix */
15358965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
15368965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
15378965ea79SLois Curfman McInnes   }
1538f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1539f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1540be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1541878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
15428965ea79SLois Curfman McInnes   A = *newmat;
15438965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
15448965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
15458965ea79SLois Curfman McInnes   }
15468965ea79SLois Curfman McInnes 
15478965ea79SLois Curfman McInnes   if (!rank) {
154887828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
15498965ea79SLois Curfman McInnes 
15508965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
15518965ea79SLois Curfman McInnes     nz = procsnz[0];
15520752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
15538965ea79SLois Curfman McInnes 
15548965ea79SLois Curfman McInnes     /* insert into matrix */
15558965ea79SLois Curfman McInnes     jj      = rstart;
15568965ea79SLois Curfman McInnes     smycols = mycols;
15578965ea79SLois Curfman McInnes     svals   = vals;
15588965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15598965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15608965ea79SLois Curfman McInnes       smycols += ourlens[i];
15618965ea79SLois Curfman McInnes       svals   += ourlens[i];
15628965ea79SLois Curfman McInnes       jj++;
15638965ea79SLois Curfman McInnes     }
15648965ea79SLois Curfman McInnes 
15658965ea79SLois Curfman McInnes     /* read in other processors and ship out */
15668965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
15678965ea79SLois Curfman McInnes       nz   = procsnz[i];
15680752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1569ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
15708965ea79SLois Curfman McInnes     }
1571606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
15723a40ed3dSBarry Smith   } else {
15738965ea79SLois Curfman McInnes     /* receive numeric values */
157484d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
15758965ea79SLois Curfman McInnes 
15768965ea79SLois Curfman McInnes     /* receive message of values*/
1577ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1578ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
157929bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
15808965ea79SLois Curfman McInnes 
15818965ea79SLois Curfman McInnes     /* insert into matrix */
15828965ea79SLois Curfman McInnes     jj      = rstart;
15838965ea79SLois Curfman McInnes     smycols = mycols;
15848965ea79SLois Curfman McInnes     svals   = vals;
15858965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15868965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15878965ea79SLois Curfman McInnes       smycols += ourlens[i];
15888965ea79SLois Curfman McInnes       svals   += ourlens[i];
15898965ea79SLois Curfman McInnes       jj++;
15908965ea79SLois Curfman McInnes     }
15918965ea79SLois Curfman McInnes   }
1592606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1593606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1594606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1595606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
15968965ea79SLois Curfman McInnes 
15976d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15986d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15993a40ed3dSBarry Smith   PetscFunctionReturn(0);
16008965ea79SLois Curfman McInnes }
160190ace30eSBarry Smith 
160290ace30eSBarry Smith 
160390ace30eSBarry Smith 
160490ace30eSBarry Smith 
160590ace30eSBarry Smith 
1606