xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 6e4ee0c6822e18ef7a0447f687dc25762b709ce5)
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__ "MatCholeskyFactorSymbolic_MPIDense"
226017fc9fSHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_MPIDense(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
230ad19415SHong Zhang {
240ad19415SHong Zhang   PetscErrorCode ierr;
250ad19415SHong Zhang 
260ad19415SHong Zhang   PetscFunctionBegin;
276017fc9fSHong Zhang   ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
280ad19415SHong Zhang   PetscFunctionReturn(0);
290ad19415SHong Zhang }
300ad19415SHong Zhang 
310ad19415SHong Zhang #undef __FUNCT__
32ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
33ab92ecdeSBarry Smith /*@
34ab92ecdeSBarry Smith 
35ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
36ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
37ab92ecdeSBarry Smith 
38ab92ecdeSBarry Smith     Input Parameter:
39ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
40ab92ecdeSBarry Smith 
41ab92ecdeSBarry Smith     Output Parameter:
42ab92ecdeSBarry Smith .      B - the inner matrix
43ab92ecdeSBarry Smith 
448e6c10adSSatish Balay     Level: intermediate
458e6c10adSSatish Balay 
46ab92ecdeSBarry Smith @*/
47ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
48ab92ecdeSBarry Smith {
49ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
50ab92ecdeSBarry Smith   PetscErrorCode ierr;
51ab92ecdeSBarry Smith   PetscTruth     flg;
52ab92ecdeSBarry Smith   PetscMPIInt    size;
53ab92ecdeSBarry Smith 
54ab92ecdeSBarry Smith   PetscFunctionBegin;
55ab92ecdeSBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
56ab92ecdeSBarry Smith   if (!flg) {  /* this check sucks! */
57ab92ecdeSBarry Smith     ierr = PetscTypeCompare((PetscObject)A,MATDENSE,&flg);CHKERRQ(ierr);
58ab92ecdeSBarry Smith     if (flg) {
59ab92ecdeSBarry Smith       ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
60ab92ecdeSBarry Smith       if (size == 1) flg = PETSC_FALSE;
61ab92ecdeSBarry Smith     }
62ab92ecdeSBarry Smith   }
63ab92ecdeSBarry Smith   if (flg) {
64ab92ecdeSBarry Smith     *B = mat->A;
65ab92ecdeSBarry Smith   } else {
66ab92ecdeSBarry Smith     *B = A;
67ab92ecdeSBarry Smith   }
68ab92ecdeSBarry Smith   PetscFunctionReturn(0);
69ab92ecdeSBarry Smith }
70ab92ecdeSBarry Smith 
71ab92ecdeSBarry Smith #undef __FUNCT__
72ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
73ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
74ba8c8a56SBarry Smith {
75ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
76ba8c8a56SBarry Smith   PetscErrorCode ierr;
77899cda47SBarry Smith   PetscInt       lrow,rstart = A->rmap.rstart,rend = A->rmap.rend;
78ba8c8a56SBarry Smith 
79ba8c8a56SBarry Smith   PetscFunctionBegin;
80ba8c8a56SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
81ba8c8a56SBarry Smith   lrow = row - rstart;
82ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
83ba8c8a56SBarry Smith   PetscFunctionReturn(0);
84ba8c8a56SBarry Smith }
85ba8c8a56SBarry Smith 
86ba8c8a56SBarry Smith #undef __FUNCT__
87ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
88ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
89ba8c8a56SBarry Smith {
90ba8c8a56SBarry Smith   PetscErrorCode ierr;
91ba8c8a56SBarry Smith 
92ba8c8a56SBarry Smith   PetscFunctionBegin;
93ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
94ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
95ba8c8a56SBarry Smith   PetscFunctionReturn(0);
96ba8c8a56SBarry Smith }
97ba8c8a56SBarry Smith 
980de54da6SSatish Balay EXTERN_C_BEGIN
994a2ae208SSatish Balay #undef __FUNCT__
1004a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
101be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
1020de54da6SSatish Balay {
1030de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
1046849ba73SBarry Smith   PetscErrorCode ierr;
105899cda47SBarry Smith   PetscInt       m = A->rmap.n,rstart = A->rmap.rstart;
10687828ca2SBarry Smith   PetscScalar    *array;
1070de54da6SSatish Balay   MPI_Comm       comm;
1080de54da6SSatish Balay 
1090de54da6SSatish Balay   PetscFunctionBegin;
110899cda47SBarry Smith   if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
1110de54da6SSatish Balay 
1120de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
1130de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
1140de54da6SSatish Balay 
1150de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
1160de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
117f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,B);CHKERRQ(ierr);
118f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr);
1195c5985e7SKris Buschelman   ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr);
1205c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
1210de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
1220de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1230de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1240de54da6SSatish Balay 
1250de54da6SSatish Balay   *iscopy = PETSC_TRUE;
1260de54da6SSatish Balay   PetscFunctionReturn(0);
1270de54da6SSatish Balay }
1280de54da6SSatish Balay EXTERN_C_END
1290de54da6SSatish Balay 
1304a2ae208SSatish Balay #undef __FUNCT__
1314a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
13213f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1338965ea79SLois Curfman McInnes {
13439b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
135dfbe8321SBarry Smith   PetscErrorCode ierr;
136899cda47SBarry Smith   PetscInt       i,j,rstart = mat->rmap.rstart,rend = mat->rmap.rend,row;
137273d9f13SBarry Smith   PetscTruth     roworiented = A->roworiented;
1388965ea79SLois Curfman McInnes 
1393a40ed3dSBarry Smith   PetscFunctionBegin;
1408965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1415ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
142899cda47SBarry Smith     if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1438965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1448965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
14539b7565bSBarry Smith       if (roworiented) {
14639b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1473a40ed3dSBarry Smith       } else {
1488965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1495ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
150899cda47SBarry Smith           if (idxn[j] >= mat->cmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
15139b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
15239b7565bSBarry Smith         }
1538965ea79SLois Curfman McInnes       }
1543a40ed3dSBarry Smith     } else {
1553782ba37SSatish Balay       if (!A->donotstash) {
15639b7565bSBarry Smith         if (roworiented) {
1578798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
158d36fbae8SSatish Balay         } else {
1598798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
16039b7565bSBarry Smith         }
161b49de8d1SLois Curfman McInnes       }
162b49de8d1SLois Curfman McInnes     }
1633782ba37SSatish Balay   }
1643a40ed3dSBarry Smith   PetscFunctionReturn(0);
165b49de8d1SLois Curfman McInnes }
166b49de8d1SLois Curfman McInnes 
1674a2ae208SSatish Balay #undef __FUNCT__
1684a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
16913f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
170b49de8d1SLois Curfman McInnes {
171b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
172dfbe8321SBarry Smith   PetscErrorCode ierr;
173899cda47SBarry Smith   PetscInt       i,j,rstart = mat->rmap.rstart,rend = mat->rmap.rend,row;
174b49de8d1SLois Curfman McInnes 
1753a40ed3dSBarry Smith   PetscFunctionBegin;
176b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
17729bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
178899cda47SBarry Smith     if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
179b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
180b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
181b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
18229bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
183899cda47SBarry Smith         if (idxn[j] >= mat->cmap.N) {
18429bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
185a8c6a408SBarry Smith         }
186b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
187b49de8d1SLois Curfman McInnes       }
188a8c6a408SBarry Smith     } else {
18929bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
1908965ea79SLois Curfman McInnes     }
1918965ea79SLois Curfman McInnes   }
1923a40ed3dSBarry Smith   PetscFunctionReturn(0);
1938965ea79SLois Curfman McInnes }
1948965ea79SLois Curfman McInnes 
1954a2ae208SSatish Balay #undef __FUNCT__
1964a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
197dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
198ff14e315SSatish Balay {
199ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
200dfbe8321SBarry Smith   PetscErrorCode ierr;
201ff14e315SSatish Balay 
2023a40ed3dSBarry Smith   PetscFunctionBegin;
203ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
2043a40ed3dSBarry Smith   PetscFunctionReturn(0);
205ff14e315SSatish Balay }
206ff14e315SSatish Balay 
2074a2ae208SSatish Balay #undef __FUNCT__
2084a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
20913f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
210ca3fa75bSLois Curfman McInnes {
211ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data,*newmatd;
212ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
2136849ba73SBarry Smith   PetscErrorCode ierr;
21413f74950SBarry Smith   PetscInt       i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
21587828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
216ca3fa75bSLois Curfman McInnes   Mat            newmat;
217ca3fa75bSLois Curfman McInnes 
218ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
219ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
220ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
221b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
222b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
223ca3fa75bSLois Curfman McInnes 
224ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
2257eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
2267eba5e9cSLois Curfman McInnes      original matrix! */
227ca3fa75bSLois Curfman McInnes 
228ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
2297eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
230ca3fa75bSLois Curfman McInnes 
231ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
232ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
23329bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2347eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
235ca3fa75bSLois Curfman McInnes     newmat = *B;
236ca3fa75bSLois Curfman McInnes   } else {
237ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
238f69a0ea3SMatthew Knepley     ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr);
239f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr);
240878740d9SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
241878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
242ca3fa75bSLois Curfman McInnes   }
243ca3fa75bSLois Curfman McInnes 
244ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
245ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
246ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
247ca3fa75bSLois Curfman McInnes 
248ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
249ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
250ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2517eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
252ca3fa75bSLois Curfman McInnes     }
253ca3fa75bSLois Curfman McInnes   }
254ca3fa75bSLois Curfman McInnes 
255ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
256ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258ca3fa75bSLois Curfman McInnes 
259ca3fa75bSLois Curfman McInnes   /* Free work space */
260ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
261ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
262ca3fa75bSLois Curfman McInnes   *B = newmat;
263ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
264ca3fa75bSLois Curfman McInnes }
265ca3fa75bSLois Curfman McInnes 
2664a2ae208SSatish Balay #undef __FUNCT__
2674a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
268dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
269ff14e315SSatish Balay {
2703a40ed3dSBarry Smith   PetscFunctionBegin;
2713a40ed3dSBarry Smith   PetscFunctionReturn(0);
272ff14e315SSatish Balay }
273ff14e315SSatish Balay 
2744a2ae208SSatish Balay #undef __FUNCT__
2754a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
276dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2778965ea79SLois Curfman McInnes {
27839ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
2798965ea79SLois Curfman McInnes   MPI_Comm       comm = mat->comm;
280dfbe8321SBarry Smith   PetscErrorCode ierr;
28113f74950SBarry Smith   PetscInt       nstash,reallocs;
2828965ea79SLois Curfman McInnes   InsertMode     addv;
2838965ea79SLois Curfman McInnes 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
2858965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
286ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
2877056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
28829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
2898965ea79SLois Curfman McInnes   }
290e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2918965ea79SLois Curfman McInnes 
292899cda47SBarry Smith   ierr = MatStashScatterBegin_Private(&mat->stash,mat->rmap.range);CHKERRQ(ierr);
2938798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
294ae15b995SBarry Smith   ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
2953a40ed3dSBarry Smith   PetscFunctionReturn(0);
2968965ea79SLois Curfman McInnes }
2978965ea79SLois Curfman McInnes 
2984a2ae208SSatish Balay #undef __FUNCT__
2994a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
300dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
3018965ea79SLois Curfman McInnes {
30239ddd567SLois Curfman McInnes   Mat_MPIDense    *mdn=(Mat_MPIDense*)mat->data;
3036849ba73SBarry Smith   PetscErrorCode  ierr;
30413f74950SBarry Smith   PetscInt        i,*row,*col,flg,j,rstart,ncols;
30513f74950SBarry Smith   PetscMPIInt     n;
30687828ca2SBarry Smith   PetscScalar     *val;
307e0fa3b82SLois Curfman McInnes   InsertMode      addv=mat->insertmode;
3088965ea79SLois Curfman McInnes 
3093a40ed3dSBarry Smith   PetscFunctionBegin;
3108965ea79SLois Curfman McInnes   /*  wait on receives */
3117ef1d9bdSSatish Balay   while (1) {
3128798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
3137ef1d9bdSSatish Balay     if (!flg) break;
3148965ea79SLois Curfman McInnes 
3157ef1d9bdSSatish Balay     for (i=0; i<n;) {
3167ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
3177ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
3187ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
3197ef1d9bdSSatish Balay       else       ncols = n-i;
3207ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
3217ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
3227ef1d9bdSSatish Balay       i = j;
3238965ea79SLois Curfman McInnes     }
3247ef1d9bdSSatish Balay   }
3258798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
3268965ea79SLois Curfman McInnes 
32739ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
32839ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3298965ea79SLois Curfman McInnes 
3306d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
33139ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3328965ea79SLois Curfman McInnes   }
3333a40ed3dSBarry Smith   PetscFunctionReturn(0);
3348965ea79SLois Curfman McInnes }
3358965ea79SLois Curfman McInnes 
3364a2ae208SSatish Balay #undef __FUNCT__
3374a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
338dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3398965ea79SLois Curfman McInnes {
340dfbe8321SBarry Smith   PetscErrorCode ierr;
34139ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3423a40ed3dSBarry Smith 
3433a40ed3dSBarry Smith   PetscFunctionBegin;
3443a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3453a40ed3dSBarry Smith   PetscFunctionReturn(0);
3468965ea79SLois Curfman McInnes }
3478965ea79SLois Curfman McInnes 
3488965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3498965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3508965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3513501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3528965ea79SLois Curfman McInnes    routine.
3538965ea79SLois Curfman McInnes */
3544a2ae208SSatish Balay #undef __FUNCT__
3554a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
356f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
3578965ea79SLois Curfman McInnes {
35839ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3596849ba73SBarry Smith   PetscErrorCode ierr;
360899cda47SBarry Smith   PetscInt       i,*owners = A->rmap.range;
36113f74950SBarry Smith   PetscInt       *nprocs,j,idx,nsends;
36213f74950SBarry Smith   PetscInt       nmax,*svalues,*starts,*owner,nrecvs;
36313f74950SBarry Smith   PetscInt       *rvalues,tag = A->tag,count,base,slen,*source;
36413f74950SBarry Smith   PetscInt       *lens,*lrows,*values;
36513f74950SBarry Smith   PetscMPIInt    n,imdex,rank = l->rank,size = l->size;
3668965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
3678965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
3688965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
36935d8aa7fSBarry Smith   PetscTruth     found;
3708965ea79SLois Curfman McInnes 
3713a40ed3dSBarry Smith   PetscFunctionBegin;
3728965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
37313f74950SBarry Smith   ierr  = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
37413f74950SBarry Smith   ierr  = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
37513f74950SBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/
3768965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3778965ea79SLois Curfman McInnes     idx = rows[i];
37835d8aa7fSBarry Smith     found = PETSC_FALSE;
3798965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3808965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
381c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3828965ea79SLois Curfman McInnes       }
3838965ea79SLois Curfman McInnes     }
38429bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3858965ea79SLois Curfman McInnes   }
386c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3878965ea79SLois Curfman McInnes 
3888965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
389c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3908965ea79SLois Curfman McInnes 
3918965ea79SLois Curfman McInnes   /* post receives:   */
39213f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
393b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3948965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
39513f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3968965ea79SLois Curfman McInnes   }
3978965ea79SLois Curfman McInnes 
3988965ea79SLois Curfman McInnes   /* do sends:
3998965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
4008965ea79SLois Curfman McInnes          the ith processor
4018965ea79SLois Curfman McInnes   */
40213f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
403b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
40413f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
4058965ea79SLois Curfman McInnes   starts[0]  = 0;
406c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
4078965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
4088965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
4098965ea79SLois Curfman McInnes   }
4108965ea79SLois Curfman McInnes 
4118965ea79SLois Curfman McInnes   starts[0] = 0;
412c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
4138965ea79SLois Curfman McInnes   count = 0;
4148965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
415c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
41613f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
4178965ea79SLois Curfman McInnes     }
4188965ea79SLois Curfman McInnes   }
419606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
4208965ea79SLois Curfman McInnes 
4218965ea79SLois Curfman McInnes   base = owners[rank];
4228965ea79SLois Curfman McInnes 
4238965ea79SLois Curfman McInnes   /*  wait on receives */
42413f74950SBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
4258965ea79SLois Curfman McInnes   source = lens + nrecvs;
4268965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
4278965ea79SLois Curfman McInnes   while (count) {
428ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4298965ea79SLois Curfman McInnes     /* unpack receives into our local space */
43013f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4318965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
4328965ea79SLois Curfman McInnes     lens[imdex]    = n;
4338965ea79SLois Curfman McInnes     slen += n;
4348965ea79SLois Curfman McInnes     count--;
4358965ea79SLois Curfman McInnes   }
436606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4378965ea79SLois Curfman McInnes 
4388965ea79SLois Curfman McInnes   /* move the data into the send scatter */
43913f74950SBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
4408965ea79SLois Curfman McInnes   count = 0;
4418965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4428965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4438965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4448965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4458965ea79SLois Curfman McInnes     }
4468965ea79SLois Curfman McInnes   }
447606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
448606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
449606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
450606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
4518965ea79SLois Curfman McInnes 
4528965ea79SLois Curfman McInnes   /* actually zap the local rows */
453f4df32b1SMatthew Knepley   ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr);
454606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4558965ea79SLois Curfman McInnes 
4568965ea79SLois Curfman McInnes   /* wait on sends */
4578965ea79SLois Curfman McInnes   if (nsends) {
458b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
459ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
460606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4618965ea79SLois Curfman McInnes   }
462606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
463606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4648965ea79SLois Curfman McInnes 
4653a40ed3dSBarry Smith   PetscFunctionReturn(0);
4668965ea79SLois Curfman McInnes }
4678965ea79SLois Curfman McInnes 
4684a2ae208SSatish Balay #undef __FUNCT__
4694a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
470dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4718965ea79SLois Curfman McInnes {
47239ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
473dfbe8321SBarry Smith   PetscErrorCode ierr;
474c456f294SBarry Smith 
4753a40ed3dSBarry Smith   PetscFunctionBegin;
47643a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
47743a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
47844cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4793a40ed3dSBarry Smith   PetscFunctionReturn(0);
4808965ea79SLois Curfman McInnes }
4818965ea79SLois Curfman McInnes 
4824a2ae208SSatish Balay #undef __FUNCT__
4834a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
484dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4858965ea79SLois Curfman McInnes {
48639ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
487dfbe8321SBarry Smith   PetscErrorCode ierr;
488c456f294SBarry Smith 
4893a40ed3dSBarry Smith   PetscFunctionBegin;
49043a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
49143a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
49244cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4933a40ed3dSBarry Smith   PetscFunctionReturn(0);
4948965ea79SLois Curfman McInnes }
4958965ea79SLois Curfman McInnes 
4964a2ae208SSatish Balay #undef __FUNCT__
4974a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
498dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
499096963f5SLois Curfman McInnes {
500096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
501dfbe8321SBarry Smith   PetscErrorCode ierr;
50287828ca2SBarry Smith   PetscScalar    zero = 0.0;
503096963f5SLois Curfman McInnes 
5043a40ed3dSBarry Smith   PetscFunctionBegin;
5052dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
5067c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
507537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
508537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
5093a40ed3dSBarry Smith   PetscFunctionReturn(0);
510096963f5SLois Curfman McInnes }
511096963f5SLois Curfman McInnes 
5124a2ae208SSatish Balay #undef __FUNCT__
5134a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
514dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
515096963f5SLois Curfman McInnes {
516096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
517dfbe8321SBarry Smith   PetscErrorCode ierr;
518096963f5SLois Curfman McInnes 
5193a40ed3dSBarry Smith   PetscFunctionBegin;
5203501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
5217c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
522537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
523537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
5243a40ed3dSBarry Smith   PetscFunctionReturn(0);
525096963f5SLois Curfman McInnes }
526096963f5SLois Curfman McInnes 
5274a2ae208SSatish Balay #undef __FUNCT__
5284a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
529dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5308965ea79SLois Curfman McInnes {
53139ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
532096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
533dfbe8321SBarry Smith   PetscErrorCode ierr;
534899cda47SBarry Smith   PetscInt       len,i,n,m = A->rmap.n,radd;
53587828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
536ed3cc1f0SBarry Smith 
5373a40ed3dSBarry Smith   PetscFunctionBegin;
5382dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5391ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
540096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
541899cda47SBarry Smith   if (n != A->rmap.N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
542899cda47SBarry Smith   len  = PetscMin(a->A->rmap.n,a->A->cmap.n);
543899cda47SBarry Smith   radd = A->rmap.rstart*m;
54444cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
545096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
546096963f5SLois Curfman McInnes   }
5471ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5483a40ed3dSBarry Smith   PetscFunctionReturn(0);
5498965ea79SLois Curfman McInnes }
5508965ea79SLois Curfman McInnes 
5514a2ae208SSatish Balay #undef __FUNCT__
5524a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
553dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5548965ea79SLois Curfman McInnes {
5553501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
556dfbe8321SBarry Smith   PetscErrorCode ierr;
557ed3cc1f0SBarry Smith 
5583a40ed3dSBarry Smith   PetscFunctionBegin;
55994d884c6SBarry Smith 
560aa482453SBarry Smith #if defined(PETSC_USE_LOG)
561899cda47SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap.N,mat->cmap.N);
5628965ea79SLois Curfman McInnes #endif
5638798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
5643501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
56505b42c5fSBarry Smith   if (mdn->lvec)   {ierr = VecDestroy(mdn->lvec);CHKERRQ(ierr);}
56605b42c5fSBarry Smith   if (mdn->Mvctx)  {ierr = VecScatterDestroy(mdn->Mvctx);CHKERRQ(ierr);}
567622d7880SLois Curfman McInnes   if (mdn->factor) {
56805b42c5fSBarry Smith     ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);
56905b42c5fSBarry Smith     ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);
57005b42c5fSBarry Smith     ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);
571606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
572622d7880SLois Curfman McInnes   }
573606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
574901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
575901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5763a40ed3dSBarry Smith   PetscFunctionReturn(0);
5778965ea79SLois Curfman McInnes }
57839ddd567SLois Curfman McInnes 
5794a2ae208SSatish Balay #undef __FUNCT__
5804a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5816849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5828965ea79SLois Curfman McInnes {
58339ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
584dfbe8321SBarry Smith   PetscErrorCode ierr;
5857056b6fcSBarry Smith 
5863a40ed3dSBarry Smith   PetscFunctionBegin;
58739ddd567SLois Curfman McInnes   if (mdn->size == 1) {
58839ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5898965ea79SLois Curfman McInnes   }
59029bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5913a40ed3dSBarry Smith   PetscFunctionReturn(0);
5928965ea79SLois Curfman McInnes }
5938965ea79SLois Curfman McInnes 
5944a2ae208SSatish Balay #undef __FUNCT__
5954a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
5966849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5978965ea79SLois Curfman McInnes {
59839ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
599dfbe8321SBarry Smith   PetscErrorCode    ierr;
60032dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
601b0a32e0cSBarry Smith   PetscViewerType   vtype;
60232077d6dSBarry Smith   PetscTruth        iascii,isdraw;
603b0a32e0cSBarry Smith   PetscViewer       sviewer;
604f3ef73ceSBarry Smith   PetscViewerFormat format;
6058965ea79SLois Curfman McInnes 
6063a40ed3dSBarry Smith   PetscFunctionBegin;
60732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
608fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
60932077d6dSBarry Smith   if (iascii) {
610b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
611b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
612456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6134e220ebcSLois Curfman McInnes       MatInfo info;
614888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
615899cda47SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap.n,
61677431f27SBarry Smith                    (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr);
617b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6183501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
6193a40ed3dSBarry Smith       PetscFunctionReturn(0);
620fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6213a40ed3dSBarry Smith       PetscFunctionReturn(0);
6228965ea79SLois Curfman McInnes     }
623f1af5d2fSBarry Smith   } else if (isdraw) {
624b0a32e0cSBarry Smith     PetscDraw  draw;
625f1af5d2fSBarry Smith     PetscTruth isnull;
626f1af5d2fSBarry Smith 
627b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
628b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
629f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
630f1af5d2fSBarry Smith   }
63177ed5343SBarry Smith 
6328965ea79SLois Curfman McInnes   if (size == 1) {
63339ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6343a40ed3dSBarry Smith   } else {
6358965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
6368965ea79SLois Curfman McInnes     Mat         A;
637899cda47SBarry Smith     PetscInt    M = mat->rmap.N,N = mat->cmap.N,m,row,i,nz;
638ba8c8a56SBarry Smith     PetscInt    *cols;
639ba8c8a56SBarry Smith     PetscScalar *vals;
6408965ea79SLois Curfman McInnes 
641f69a0ea3SMatthew Knepley     ierr = MatCreate(mat->comm,&A);CHKERRQ(ierr);
6428965ea79SLois Curfman McInnes     if (!rank) {
643f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
6443a40ed3dSBarry Smith     } else {
645f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
6468965ea79SLois Curfman McInnes     }
647be5d1d56SKris Buschelman     /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */
648878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
649878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
65052e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
6518965ea79SLois Curfman McInnes 
65239ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
65339ddd567SLois Curfman McInnes        but it's quick for now */
65451022da4SBarry Smith     A->insertmode = INSERT_VALUES;
655899cda47SBarry Smith     row = mat->rmap.rstart; m = mdn->A->rmap.n;
6568965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
657ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
658ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
659ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
66039ddd567SLois Curfman McInnes       row++;
6618965ea79SLois Curfman McInnes     }
6628965ea79SLois Curfman McInnes 
6636d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6646d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
665b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
666b9b97703SBarry Smith     if (!rank) {
6676831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
6688965ea79SLois Curfman McInnes     }
669b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
670b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6718965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
6728965ea79SLois Curfman McInnes   }
6733a40ed3dSBarry Smith   PetscFunctionReturn(0);
6748965ea79SLois Curfman McInnes }
6758965ea79SLois Curfman McInnes 
6764a2ae208SSatish Balay #undef __FUNCT__
6774a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
678dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
6798965ea79SLois Curfman McInnes {
680dfbe8321SBarry Smith   PetscErrorCode ierr;
68132077d6dSBarry Smith   PetscTruth     iascii,isbinary,isdraw,issocket;
6828965ea79SLois Curfman McInnes 
683433994e6SBarry Smith   PetscFunctionBegin;
6840f5bd95cSBarry Smith 
68532077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
686fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
687b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
688fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6890f5bd95cSBarry Smith 
69032077d6dSBarry Smith   if (iascii || issocket || isdraw) {
691f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
6920f5bd95cSBarry Smith   } else if (isbinary) {
6933a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6945cd90555SBarry Smith   } else {
695958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6968965ea79SLois Curfman McInnes   }
6973a40ed3dSBarry Smith   PetscFunctionReturn(0);
6988965ea79SLois Curfman McInnes }
6998965ea79SLois Curfman McInnes 
7004a2ae208SSatish Balay #undef __FUNCT__
7014a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
702dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7038965ea79SLois Curfman McInnes {
7043501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7053501a2bdSLois Curfman McInnes   Mat            mdn = mat->A;
706dfbe8321SBarry Smith   PetscErrorCode ierr;
707329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7088965ea79SLois Curfman McInnes 
7093a40ed3dSBarry Smith   PetscFunctionBegin;
710899cda47SBarry Smith   info->rows_global    = (double)A->rmap.N;
711899cda47SBarry Smith   info->columns_global = (double)A->cmap.N;
712899cda47SBarry Smith   info->rows_local     = (double)A->rmap.n;
713899cda47SBarry Smith   info->columns_local  = (double)A->cmap.N;
7144e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
7154e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7164e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7174e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7188965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7194e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7204e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7214e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7224e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7234e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7248965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
725d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
7264e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7274e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7284e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7294e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7304e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7318965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
732d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
7334e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7344e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7354e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7364e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7374e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7388965ea79SLois Curfman McInnes   }
7394e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
7404e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
7414e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
7423a40ed3dSBarry Smith   PetscFunctionReturn(0);
7438965ea79SLois Curfman McInnes }
7448965ea79SLois Curfman McInnes 
7454a2ae208SSatish Balay #undef __FUNCT__
7464a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
747dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op)
7488965ea79SLois Curfman McInnes {
74939ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
750dfbe8321SBarry Smith   PetscErrorCode ierr;
7518965ea79SLois Curfman McInnes 
7523a40ed3dSBarry Smith   PetscFunctionBegin;
75312c028f9SKris Buschelman   switch (op) {
75412c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
75512c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
75612c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
75712c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
75812c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
75912c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
760273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
76112c028f9SKris Buschelman     break;
76212c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
763273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
764273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
76512c028f9SKris Buschelman     break;
76612c028f9SKris Buschelman   case MAT_ROWS_SORTED:
76712c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
76812c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
76912c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
770ad86a440SBarry Smith     ierr = PetscInfo1(A,"Option %d ignored\n",op);CHKERRQ(ierr);
77112c028f9SKris Buschelman     break;
77212c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
773273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
774273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
77512c028f9SKris Buschelman     break;
77612c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
777273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
77812c028f9SKris Buschelman     break;
77912c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
78029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
78177e54ba9SKris Buschelman   case MAT_SYMMETRIC:
78277e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
7839a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
7849a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
7859a4540c5SBarry Smith   case MAT_HERMITIAN:
7869a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
7879a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
7889a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
78977e54ba9SKris Buschelman     break;
79012c028f9SKris Buschelman   default:
791ad86a440SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
7923a40ed3dSBarry Smith   }
7933a40ed3dSBarry Smith   PetscFunctionReturn(0);
7948965ea79SLois Curfman McInnes }
7958965ea79SLois Curfman McInnes 
7968965ea79SLois Curfman McInnes 
7974a2ae208SSatish Balay #undef __FUNCT__
7984a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
799dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8005b2fa520SLois Curfman McInnes {
8015b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8025b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
80387828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
804dfbe8321SBarry Smith   PetscErrorCode ierr;
805899cda47SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap.n,n=mdn->A->cmap.n;
8065b2fa520SLois Curfman McInnes 
8075b2fa520SLois Curfman McInnes   PetscFunctionBegin;
80872d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8095b2fa520SLois Curfman McInnes   if (ll) {
81072d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
811175be7b4SMatthew Knepley     if (s2a != s2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8121ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8135b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8145b2fa520SLois Curfman McInnes       x = l[i];
8155b2fa520SLois Curfman McInnes       v = mat->v + i;
8165b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8175b2fa520SLois Curfman McInnes     }
8181ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
819efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8205b2fa520SLois Curfman McInnes   }
8215b2fa520SLois Curfman McInnes   if (rr) {
822175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
823175be7b4SMatthew Knepley     if (s3a != s3) SETERRQ2(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
8245b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8255b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
8261ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8275b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8285b2fa520SLois Curfman McInnes       x = r[i];
8295b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8305b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
8315b2fa520SLois Curfman McInnes     }
8321ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
833efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8345b2fa520SLois Curfman McInnes   }
8355b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8365b2fa520SLois Curfman McInnes }
8375b2fa520SLois Curfman McInnes 
8384a2ae208SSatish Balay #undef __FUNCT__
8394a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
840dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
841096963f5SLois Curfman McInnes {
8423501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8433501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
844dfbe8321SBarry Smith   PetscErrorCode ierr;
84513f74950SBarry Smith   PetscInt       i,j;
846329f5518SBarry Smith   PetscReal      sum = 0.0;
84787828ca2SBarry Smith   PetscScalar    *v = mat->v;
8483501a2bdSLois Curfman McInnes 
8493a40ed3dSBarry Smith   PetscFunctionBegin;
8503501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
851064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
8523501a2bdSLois Curfman McInnes   } else {
8533501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
854899cda47SBarry Smith       for (i=0; i<mdn->A->cmap.n*mdn->A->rmap.n; i++) {
855aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
856329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
8573501a2bdSLois Curfman McInnes #else
8583501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
8593501a2bdSLois Curfman McInnes #endif
8603501a2bdSLois Curfman McInnes       }
861064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
862064f8208SBarry Smith       *nrm = sqrt(*nrm);
863899cda47SBarry Smith       ierr = PetscLogFlops(2*mdn->A->cmap.n*mdn->A->rmap.n);CHKERRQ(ierr);
8643a40ed3dSBarry Smith     } else if (type == NORM_1) {
865329f5518SBarry Smith       PetscReal *tmp,*tmp2;
866899cda47SBarry Smith       ierr = PetscMalloc(2*A->cmap.N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
867899cda47SBarry Smith       tmp2 = tmp + A->cmap.N;
868899cda47SBarry Smith       ierr = PetscMemzero(tmp,2*A->cmap.N*sizeof(PetscReal));CHKERRQ(ierr);
869064f8208SBarry Smith       *nrm = 0.0;
8703501a2bdSLois Curfman McInnes       v = mat->v;
871899cda47SBarry Smith       for (j=0; j<mdn->A->cmap.n; j++) {
872899cda47SBarry Smith         for (i=0; i<mdn->A->rmap.n; i++) {
87367e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8743501a2bdSLois Curfman McInnes         }
8753501a2bdSLois Curfman McInnes       }
876899cda47SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->cmap.N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
877899cda47SBarry Smith       for (j=0; j<A->cmap.N; j++) {
878064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8793501a2bdSLois Curfman McInnes       }
880606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
881899cda47SBarry Smith       ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
8823a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
883329f5518SBarry Smith       PetscReal ntemp;
8843501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
885064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
8863a40ed3dSBarry Smith     } else {
88729bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8883501a2bdSLois Curfman McInnes     }
8893501a2bdSLois Curfman McInnes   }
8903a40ed3dSBarry Smith   PetscFunctionReturn(0);
8913501a2bdSLois Curfman McInnes }
8923501a2bdSLois Curfman McInnes 
8934a2ae208SSatish Balay #undef __FUNCT__
8944a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
895dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout)
8963501a2bdSLois Curfman McInnes {
8973501a2bdSLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
8983501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
8993501a2bdSLois Curfman McInnes   Mat            B;
900899cda47SBarry Smith   PetscInt       M = A->rmap.N,N = A->cmap.N,m,n,*rwork,rstart = A->rmap.rstart;
9016849ba73SBarry Smith   PetscErrorCode ierr;
90213f74950SBarry Smith   PetscInt       j,i;
90387828ca2SBarry Smith   PetscScalar    *v;
9043501a2bdSLois Curfman McInnes 
9053a40ed3dSBarry Smith   PetscFunctionBegin;
9067c922b88SBarry Smith   if (!matout && M != N) {
90729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
9087056b6fcSBarry Smith   }
909f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
910f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr);
911878740d9SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
912878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
9133501a2bdSLois Curfman McInnes 
914899cda47SBarry Smith   m = a->A->rmap.n; n = a->A->cmap.n; v = Aloc->v;
91513f74950SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
9163501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
9173501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
9183501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9193501a2bdSLois Curfman McInnes     v   += m;
9203501a2bdSLois Curfman McInnes   }
921606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9226d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9236d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9247c922b88SBarry Smith   if (matout) {
9253501a2bdSLois Curfman McInnes     *matout = B;
9263501a2bdSLois Curfman McInnes   } else {
927273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
9283501a2bdSLois Curfman McInnes   }
9293a40ed3dSBarry Smith   PetscFunctionReturn(0);
930096963f5SLois Curfman McInnes }
931096963f5SLois Curfman McInnes 
932d9eff348SSatish Balay #include "petscblaslapack.h"
9334a2ae208SSatish Balay #undef __FUNCT__
9344a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
935f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha)
93644cd7ae7SLois Curfman McInnes {
93744cd7ae7SLois Curfman McInnes   Mat_MPIDense   *A = (Mat_MPIDense*)inA->data;
93844cd7ae7SLois Curfman McInnes   Mat_SeqDense   *a = (Mat_SeqDense*)A->A->data;
939f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
940899cda47SBarry Smith   PetscBLASInt   one = 1,nz = (PetscBLASInt)inA->rmap.n*inA->cmap.N;
941efee365bSSatish Balay   PetscErrorCode ierr;
94244cd7ae7SLois Curfman McInnes 
9433a40ed3dSBarry Smith   PetscFunctionBegin;
944f4df32b1SMatthew Knepley   BLASscal_(&nz,&oalpha,a->v,&one);
945efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
9463a40ed3dSBarry Smith   PetscFunctionReturn(0);
94744cd7ae7SLois Curfman McInnes }
94844cd7ae7SLois Curfman McInnes 
9496849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
9508965ea79SLois Curfman McInnes 
9514a2ae208SSatish Balay #undef __FUNCT__
9524a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
953dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
954273d9f13SBarry Smith {
955dfbe8321SBarry Smith   PetscErrorCode ierr;
956273d9f13SBarry Smith 
957273d9f13SBarry Smith   PetscFunctionBegin;
958273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
959273d9f13SBarry Smith   PetscFunctionReturn(0);
960273d9f13SBarry Smith }
961273d9f13SBarry Smith 
9628965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
96309dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
96409dc0095SBarry Smith        MatGetRow_MPIDense,
96509dc0095SBarry Smith        MatRestoreRow_MPIDense,
96609dc0095SBarry Smith        MatMult_MPIDense,
96797304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
9687c922b88SBarry Smith        MatMultTranspose_MPIDense,
9697c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
9708965ea79SLois Curfman McInnes        0,
97109dc0095SBarry Smith        0,
97209dc0095SBarry Smith        0,
97397304618SKris Buschelman /*10*/ 0,
97409dc0095SBarry Smith        0,
97509dc0095SBarry Smith        0,
97609dc0095SBarry Smith        0,
97709dc0095SBarry Smith        MatTranspose_MPIDense,
97897304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
979*6e4ee0c6SHong Zhang        MatEqual_MPIDense,
98009dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9815b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
98209dc0095SBarry Smith        MatNorm_MPIDense,
98397304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
98409dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
98509dc0095SBarry Smith        0,
98609dc0095SBarry Smith        MatSetOption_MPIDense,
98709dc0095SBarry Smith        MatZeroEntries_MPIDense,
98897304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
9890ad19415SHong Zhang        MatLUFactorSymbolic_MPIDense,
990919b68f7SBarry Smith        0,
9910ad19415SHong Zhang        MatCholeskyFactorSymbolic_MPIDense,
992919b68f7SBarry Smith        0,
99397304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
994273d9f13SBarry Smith        0,
99509dc0095SBarry Smith        0,
99609dc0095SBarry Smith        MatGetArray_MPIDense,
99709dc0095SBarry Smith        MatRestoreArray_MPIDense,
99897304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
99909dc0095SBarry Smith        0,
100009dc0095SBarry Smith        0,
100109dc0095SBarry Smith        0,
100209dc0095SBarry Smith        0,
100397304618SKris Buschelman /*40*/ 0,
10042ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
100509dc0095SBarry Smith        0,
100609dc0095SBarry Smith        MatGetValues_MPIDense,
100709dc0095SBarry Smith        0,
100897304618SKris Buschelman /*45*/ 0,
100909dc0095SBarry Smith        MatScale_MPIDense,
101009dc0095SBarry Smith        0,
101109dc0095SBarry Smith        0,
101209dc0095SBarry Smith        0,
1013521d7252SBarry Smith /*50*/ 0,
101409dc0095SBarry Smith        0,
101509dc0095SBarry Smith        0,
101609dc0095SBarry Smith        0,
101709dc0095SBarry Smith        0,
101897304618SKris Buschelman /*55*/ 0,
101909dc0095SBarry Smith        0,
102009dc0095SBarry Smith        0,
102109dc0095SBarry Smith        0,
102209dc0095SBarry Smith        0,
102397304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
1024b9b97703SBarry Smith        MatDestroy_MPIDense,
1025b9b97703SBarry Smith        MatView_MPIDense,
1026357abbc8SBarry Smith        0,
102797304618SKris Buschelman        0,
102897304618SKris Buschelman /*65*/ 0,
102997304618SKris Buschelman        0,
103097304618SKris Buschelman        0,
103197304618SKris Buschelman        0,
103297304618SKris Buschelman        0,
103397304618SKris Buschelman /*70*/ 0,
103497304618SKris Buschelman        0,
103597304618SKris Buschelman        0,
103697304618SKris Buschelman        0,
103797304618SKris Buschelman        0,
103897304618SKris Buschelman /*75*/ 0,
103997304618SKris Buschelman        0,
104097304618SKris Buschelman        0,
104197304618SKris Buschelman        0,
104297304618SKris Buschelman        0,
104397304618SKris Buschelman /*80*/ 0,
104497304618SKris Buschelman        0,
104597304618SKris Buschelman        0,
104697304618SKris Buschelman        0,
1047865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense,
1048865e5f61SKris Buschelman        0,
1049865e5f61SKris Buschelman        0,
1050865e5f61SKris Buschelman        0,
1051865e5f61SKris Buschelman        0,
1052865e5f61SKris Buschelman        0,
1053865e5f61SKris Buschelman /*90*/ 0,
1054865e5f61SKris Buschelman        0,
1055865e5f61SKris Buschelman        0,
1056865e5f61SKris Buschelman        0,
1057865e5f61SKris Buschelman        0,
1058865e5f61SKris Buschelman /*95*/ 0,
1059865e5f61SKris Buschelman        0,
1060865e5f61SKris Buschelman        0,
1061865e5f61SKris Buschelman        0};
10628965ea79SLois Curfman McInnes 
1063273d9f13SBarry Smith EXTERN_C_BEGIN
10644a2ae208SSatish Balay #undef __FUNCT__
1065a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
1066be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1067a23d5eceSKris Buschelman {
1068a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1069dfbe8321SBarry Smith   PetscErrorCode ierr;
1070a23d5eceSKris Buschelman 
1071a23d5eceSKris Buschelman   PetscFunctionBegin;
1072a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1073a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1074a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1075a23d5eceSKris Buschelman 
1076a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
1077f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1078899cda47SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap.n,mat->cmap.N,mat->rmap.n,mat->cmap.N);CHKERRQ(ierr);
10795c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
10805c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
108152e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1082a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1083a23d5eceSKris Buschelman }
1084a23d5eceSKris Buschelman EXTERN_C_END
1085a23d5eceSKris Buschelman 
10860bad9183SKris Buschelman /*MC
1087fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
10880bad9183SKris Buschelman 
10890bad9183SKris Buschelman    Options Database Keys:
10900bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
10910bad9183SKris Buschelman 
10920bad9183SKris Buschelman   Level: beginner
10930bad9183SKris Buschelman 
10940bad9183SKris Buschelman .seealso: MatCreateMPIDense
10950bad9183SKris Buschelman M*/
10960bad9183SKris Buschelman 
1097a23d5eceSKris Buschelman EXTERN_C_BEGIN
1098a23d5eceSKris Buschelman #undef __FUNCT__
10994a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1100be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat)
1101273d9f13SBarry Smith {
1102273d9f13SBarry Smith   Mat_MPIDense   *a;
1103dfbe8321SBarry Smith   PetscErrorCode ierr;
1104273d9f13SBarry Smith 
1105273d9f13SBarry Smith   PetscFunctionBegin;
1106b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1107b0a32e0cSBarry Smith   mat->data         = (void*)a;
1108273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1109273d9f13SBarry Smith   mat->factor       = 0;
1110273d9f13SBarry Smith   mat->mapping      = 0;
1111273d9f13SBarry Smith 
1112273d9f13SBarry Smith   a->factor       = 0;
1113273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1114273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
1115273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
1116273d9f13SBarry Smith 
1117899cda47SBarry Smith   mat->rmap.bs = mat->cmap.bs = 1;
1118899cda47SBarry Smith   ierr = PetscMapInitialize(mat->comm,&mat->rmap);CHKERRQ(ierr);
1119899cda47SBarry Smith   ierr = PetscMapInitialize(mat->comm,&mat->cmap);CHKERRQ(ierr);
1120899cda47SBarry Smith   a->nvec = mat->cmap.n;
1121273d9f13SBarry Smith 
1122273d9f13SBarry Smith   /* build cache for off array entries formed */
1123273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1124273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1125273d9f13SBarry Smith 
1126273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1127273d9f13SBarry Smith   a->lvec        = 0;
1128273d9f13SBarry Smith   a->Mvctx       = 0;
1129273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1130273d9f13SBarry Smith 
1131273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1132273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1133273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1134a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1135a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1136a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1137273d9f13SBarry Smith   PetscFunctionReturn(0);
1138273d9f13SBarry Smith }
1139273d9f13SBarry Smith EXTERN_C_END
1140273d9f13SBarry Smith 
1141209238afSKris Buschelman /*MC
1142002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1143209238afSKris Buschelman 
1144209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1145209238afSKris Buschelman    and MATMPIDENSE otherwise.
1146209238afSKris Buschelman 
1147209238afSKris Buschelman    Options Database Keys:
1148209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1149209238afSKris Buschelman 
1150209238afSKris Buschelman   Level: beginner
1151209238afSKris Buschelman 
1152209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1153209238afSKris Buschelman M*/
1154209238afSKris Buschelman 
1155209238afSKris Buschelman EXTERN_C_BEGIN
1156209238afSKris Buschelman #undef __FUNCT__
1157209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1158be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A)
1159dfbe8321SBarry Smith {
11606849ba73SBarry Smith   PetscErrorCode ierr;
116113f74950SBarry Smith   PetscMPIInt    size;
1162209238afSKris Buschelman 
1163209238afSKris Buschelman   PetscFunctionBegin;
1164209238afSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr);
1165209238afSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1166209238afSKris Buschelman   if (size == 1) {
1167209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1168209238afSKris Buschelman   } else {
1169209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1170209238afSKris Buschelman   }
1171209238afSKris Buschelman   PetscFunctionReturn(0);
1172209238afSKris Buschelman }
1173209238afSKris Buschelman EXTERN_C_END
1174209238afSKris Buschelman 
11754a2ae208SSatish Balay #undef __FUNCT__
11764a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1177273d9f13SBarry Smith /*@C
1178273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1179273d9f13SBarry Smith 
1180273d9f13SBarry Smith    Not collective
1181273d9f13SBarry Smith 
1182273d9f13SBarry Smith    Input Parameters:
1183273d9f13SBarry Smith .  A - the matrix
1184273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1185273d9f13SBarry Smith    to control all matrix memory allocation.
1186273d9f13SBarry Smith 
1187273d9f13SBarry Smith    Notes:
1188273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1189273d9f13SBarry Smith    storage by columns.
1190273d9f13SBarry Smith 
1191273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1192273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1193273d9f13SBarry Smith    set data=PETSC_NULL.
1194273d9f13SBarry Smith 
1195273d9f13SBarry Smith    Level: intermediate
1196273d9f13SBarry Smith 
1197273d9f13SBarry Smith .keywords: matrix,dense, parallel
1198273d9f13SBarry Smith 
1199273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1200273d9f13SBarry Smith @*/
1201be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1202273d9f13SBarry Smith {
1203dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar *);
1204273d9f13SBarry Smith 
1205273d9f13SBarry Smith   PetscFunctionBegin;
1206565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1207a23d5eceSKris Buschelman   if (f) {
1208a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1209a23d5eceSKris Buschelman   }
1210273d9f13SBarry Smith   PetscFunctionReturn(0);
1211273d9f13SBarry Smith }
1212273d9f13SBarry Smith 
12134a2ae208SSatish Balay #undef __FUNCT__
12144a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
12158965ea79SLois Curfman McInnes /*@C
121639ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
12178965ea79SLois Curfman McInnes 
1218db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1219db81eaa0SLois Curfman McInnes 
12208965ea79SLois Curfman McInnes    Input Parameters:
1221db81eaa0SLois Curfman McInnes +  comm - MPI communicator
12228965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1223db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
12248965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1225db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
12267f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1227dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
12288965ea79SLois Curfman McInnes 
12298965ea79SLois Curfman McInnes    Output Parameter:
1230477f1c0bSLois Curfman McInnes .  A - the matrix
12318965ea79SLois Curfman McInnes 
1232b259b22eSLois Curfman McInnes    Notes:
123339ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
123439ddd567SLois Curfman McInnes    storage by columns.
12358965ea79SLois Curfman McInnes 
123618f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
123718f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
12387f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
123918f449edSLois Curfman McInnes 
12408965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
12418965ea79SLois Curfman McInnes    (possibly both).
12428965ea79SLois Curfman McInnes 
1243027ccd11SLois Curfman McInnes    Level: intermediate
1244027ccd11SLois Curfman McInnes 
124539ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
12468965ea79SLois Curfman McInnes 
124739ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
12488965ea79SLois Curfman McInnes @*/
1249be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
12508965ea79SLois Curfman McInnes {
12516849ba73SBarry Smith   PetscErrorCode ierr;
125213f74950SBarry Smith   PetscMPIInt    size;
12538965ea79SLois Curfman McInnes 
12543a40ed3dSBarry Smith   PetscFunctionBegin;
1255f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1256f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1257273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1258273d9f13SBarry Smith   if (size > 1) {
1259273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1260273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1261273d9f13SBarry Smith   } else {
1262273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1263273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
12648c469469SLois Curfman McInnes   }
12653a40ed3dSBarry Smith   PetscFunctionReturn(0);
12668965ea79SLois Curfman McInnes }
12678965ea79SLois Curfman McInnes 
12684a2ae208SSatish Balay #undef __FUNCT__
12694a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
12706849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
12718965ea79SLois Curfman McInnes {
12728965ea79SLois Curfman McInnes   Mat            mat;
12733501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1274dfbe8321SBarry Smith   PetscErrorCode ierr;
12758965ea79SLois Curfman McInnes 
12763a40ed3dSBarry Smith   PetscFunctionBegin;
12778965ea79SLois Curfman McInnes   *newmat       = 0;
1278f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&mat);CHKERRQ(ierr);
1279899cda47SBarry Smith   ierr = MatSetSizes(mat,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
1280be5d1d56SKris Buschelman   ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr);
1281834f8fabSBarry Smith   a                 = (Mat_MPIDense*)mat->data;
1282e04c1aa4SHong Zhang   ierr              = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
12833501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1284c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1285273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
12868965ea79SLois Curfman McInnes 
1287899cda47SBarry Smith   mat->rmap.rstart     = A->rmap.rstart;
1288899cda47SBarry Smith   mat->rmap.rend       = A->rmap.rend;
12898965ea79SLois Curfman McInnes   a->size              = oldmat->size;
12908965ea79SLois Curfman McInnes   a->rank              = oldmat->rank;
1291e0fa3b82SLois Curfman McInnes   mat->insertmode      = NOT_SET_VALUES;
1292b9b97703SBarry Smith   a->nvec              = oldmat->nvec;
12933782ba37SSatish Balay   a->donotstash        = oldmat->donotstash;
1294e04c1aa4SHong Zhang 
1295899cda47SBarry Smith   ierr = PetscMalloc((a->size+1)*sizeof(PetscInt),&mat->rmap.range);CHKERRQ(ierr);
1296899cda47SBarry Smith   ierr = PetscMalloc((a->size+1)*sizeof(PetscInt),&mat->cmap.range);CHKERRQ(ierr);
1297899cda47SBarry Smith   ierr = PetscMemcpy(mat->rmap.range,A->rmap.range,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
1298899cda47SBarry Smith   ierr = PetscMemcpy(mat->cmap.range,A->cmap.range,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
12998798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
13008965ea79SLois Curfman McInnes 
1301329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
13025609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
130352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
13048965ea79SLois Curfman McInnes   *newmat = mat;
13053a40ed3dSBarry Smith   PetscFunctionReturn(0);
13068965ea79SLois Curfman McInnes }
13078965ea79SLois Curfman McInnes 
1308e090d566SSatish Balay #include "petscsys.h"
13098965ea79SLois Curfman McInnes 
13104a2ae208SSatish Balay #undef __FUNCT__
13114a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
1312f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat)
131390ace30eSBarry Smith {
13146849ba73SBarry Smith   PetscErrorCode ierr;
131532dcc486SBarry Smith   PetscMPIInt    rank,size;
131613f74950SBarry Smith   PetscInt       *rowners,i,m,nz,j;
131787828ca2SBarry Smith   PetscScalar    *array,*vals,*vals_ptr;
131890ace30eSBarry Smith   MPI_Status     status;
131990ace30eSBarry Smith 
13203a40ed3dSBarry Smith   PetscFunctionBegin;
1321d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1322d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
132390ace30eSBarry Smith 
132490ace30eSBarry Smith   /* determine ownership of all rows */
132590ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
132613f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
132713f74950SBarry Smith   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
132890ace30eSBarry Smith   rowners[0] = 0;
132990ace30eSBarry Smith   for (i=2; i<=size; i++) {
133090ace30eSBarry Smith     rowners[i] += rowners[i-1];
133190ace30eSBarry Smith   }
133290ace30eSBarry Smith 
1333f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1334f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1335be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1336878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
133790ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
133890ace30eSBarry Smith 
133990ace30eSBarry Smith   if (!rank) {
134087828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
134190ace30eSBarry Smith 
134290ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
13430752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
134490ace30eSBarry Smith 
134590ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
134690ace30eSBarry Smith     vals_ptr = vals;
134790ace30eSBarry Smith     for (i=0; i<m; i++) {
134890ace30eSBarry Smith       for (j=0; j<N; j++) {
134990ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
135090ace30eSBarry Smith       }
135190ace30eSBarry Smith     }
135290ace30eSBarry Smith 
135390ace30eSBarry Smith     /* read in other processors and ship out */
135490ace30eSBarry Smith     for (i=1; i<size; i++) {
135590ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
13560752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1357ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
135890ace30eSBarry Smith     }
13593a40ed3dSBarry Smith   } else {
136090ace30eSBarry Smith     /* receive numeric values */
136187828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
136290ace30eSBarry Smith 
136390ace30eSBarry Smith     /* receive message of values*/
1364ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
136590ace30eSBarry Smith 
136690ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
136790ace30eSBarry Smith     vals_ptr = vals;
136890ace30eSBarry Smith     for (i=0; i<m; i++) {
136990ace30eSBarry Smith       for (j=0; j<N; j++) {
137090ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
137190ace30eSBarry Smith       }
137290ace30eSBarry Smith     }
137390ace30eSBarry Smith   }
1374606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1375606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
13766d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13776d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13783a40ed3dSBarry Smith   PetscFunctionReturn(0);
137990ace30eSBarry Smith }
138090ace30eSBarry Smith 
13814a2ae208SSatish Balay #undef __FUNCT__
13824a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
1383f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat)
13848965ea79SLois Curfman McInnes {
13858965ea79SLois Curfman McInnes   Mat            A;
138687828ca2SBarry Smith   PetscScalar    *vals,*svals;
138719bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
13888965ea79SLois Curfman McInnes   MPI_Status     status;
138913f74950SBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
139013f74950SBarry Smith   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
139113f74950SBarry Smith   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
139213f74950SBarry Smith   PetscInt       i,nz,j,rstart,rend;
139313f74950SBarry Smith   int            fd;
13946849ba73SBarry Smith   PetscErrorCode ierr;
13958965ea79SLois Curfman McInnes 
13963a40ed3dSBarry Smith   PetscFunctionBegin;
1397d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1398d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13998965ea79SLois Curfman McInnes   if (!rank) {
1400b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
14010752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1402552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
14038965ea79SLois Curfman McInnes   }
14048965ea79SLois Curfman McInnes 
140513f74950SBarry Smith   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
140690ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
140790ace30eSBarry Smith 
140890ace30eSBarry Smith   /*
140990ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
141090ace30eSBarry Smith   */
141190ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1412be5d1d56SKris Buschelman     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr);
14133a40ed3dSBarry Smith     PetscFunctionReturn(0);
141490ace30eSBarry Smith   }
141590ace30eSBarry Smith 
14168965ea79SLois Curfman McInnes   /* determine ownership of all rows */
14178965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
141813f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
1419ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
14208965ea79SLois Curfman McInnes   rowners[0] = 0;
14218965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
14228965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
14238965ea79SLois Curfman McInnes   }
14248965ea79SLois Curfman McInnes   rstart = rowners[rank];
14258965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
14268965ea79SLois Curfman McInnes 
14278965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
142813f74950SBarry Smith   ierr    = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr);
14298965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
14308965ea79SLois Curfman McInnes   if (!rank) {
143113f74950SBarry Smith     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
14320752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
143313f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
14348965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
14351466f1e1SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
1436606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1437ca161407SBarry Smith   } else {
14381466f1e1SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
14398965ea79SLois Curfman McInnes   }
14408965ea79SLois Curfman McInnes 
14418965ea79SLois Curfman McInnes   if (!rank) {
14428965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
144313f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
144413f74950SBarry Smith     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
14458965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14468965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
14478965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
14488965ea79SLois Curfman McInnes       }
14498965ea79SLois Curfman McInnes     }
1450606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
14518965ea79SLois Curfman McInnes 
14528965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
14538965ea79SLois Curfman McInnes     maxnz = 0;
14548965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
14550452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
14568965ea79SLois Curfman McInnes     }
145713f74950SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
14588965ea79SLois Curfman McInnes 
14598965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
14608965ea79SLois Curfman McInnes     nz = procsnz[0];
146113f74950SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
14620752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
14638965ea79SLois Curfman McInnes 
14648965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
14658965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14668965ea79SLois Curfman McInnes       nz   = procsnz[i];
14670752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
146813f74950SBarry Smith       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
14698965ea79SLois Curfman McInnes     }
1470606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
14713a40ed3dSBarry Smith   } else {
14728965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
14738965ea79SLois Curfman McInnes     nz = 0;
14748965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14758965ea79SLois Curfman McInnes       nz += ourlens[i];
14768965ea79SLois Curfman McInnes     }
147713f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
14788965ea79SLois Curfman McInnes 
14798965ea79SLois Curfman McInnes     /* receive message of column indices*/
148013f74950SBarry Smith     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
148113f74950SBarry Smith     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
148229bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14838965ea79SLois Curfman McInnes   }
14848965ea79SLois Curfman McInnes 
14858965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
148613f74950SBarry Smith   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
14878965ea79SLois Curfman McInnes   jj = 0;
14888965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14898965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
14908965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
14918965ea79SLois Curfman McInnes       jj++;
14928965ea79SLois Curfman McInnes     }
14938965ea79SLois Curfman McInnes   }
14948965ea79SLois Curfman McInnes 
14958965ea79SLois Curfman McInnes   /* create our matrix */
14968965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14978965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
14988965ea79SLois Curfman McInnes   }
1499f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1500f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1501be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1502878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
15038965ea79SLois Curfman McInnes   A = *newmat;
15048965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
15058965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
15068965ea79SLois Curfman McInnes   }
15078965ea79SLois Curfman McInnes 
15088965ea79SLois Curfman McInnes   if (!rank) {
150987828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
15108965ea79SLois Curfman McInnes 
15118965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
15128965ea79SLois Curfman McInnes     nz = procsnz[0];
15130752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
15148965ea79SLois Curfman McInnes 
15158965ea79SLois Curfman McInnes     /* insert into matrix */
15168965ea79SLois Curfman McInnes     jj      = rstart;
15178965ea79SLois Curfman McInnes     smycols = mycols;
15188965ea79SLois Curfman McInnes     svals   = vals;
15198965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15208965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15218965ea79SLois Curfman McInnes       smycols += ourlens[i];
15228965ea79SLois Curfman McInnes       svals   += ourlens[i];
15238965ea79SLois Curfman McInnes       jj++;
15248965ea79SLois Curfman McInnes     }
15258965ea79SLois Curfman McInnes 
15268965ea79SLois Curfman McInnes     /* read in other processors and ship out */
15278965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
15288965ea79SLois Curfman McInnes       nz   = procsnz[i];
15290752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1530ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
15318965ea79SLois Curfman McInnes     }
1532606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
15333a40ed3dSBarry Smith   } else {
15348965ea79SLois Curfman McInnes     /* receive numeric values */
153584d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
15368965ea79SLois Curfman McInnes 
15378965ea79SLois Curfman McInnes     /* receive message of values*/
1538ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1539ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
154029bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
15418965ea79SLois Curfman McInnes 
15428965ea79SLois Curfman McInnes     /* insert into matrix */
15438965ea79SLois Curfman McInnes     jj      = rstart;
15448965ea79SLois Curfman McInnes     smycols = mycols;
15458965ea79SLois Curfman McInnes     svals   = vals;
15468965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
15478965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
15488965ea79SLois Curfman McInnes       smycols += ourlens[i];
15498965ea79SLois Curfman McInnes       svals   += ourlens[i];
15508965ea79SLois Curfman McInnes       jj++;
15518965ea79SLois Curfman McInnes     }
15528965ea79SLois Curfman McInnes   }
1553606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1554606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1555606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1556606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
15578965ea79SLois Curfman McInnes 
15586d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15596d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15603a40ed3dSBarry Smith   PetscFunctionReturn(0);
15618965ea79SLois Curfman McInnes }
156290ace30eSBarry Smith 
1563*6e4ee0c6SHong Zhang #undef __FUNCT__
1564*6e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense"
1565*6e4ee0c6SHong Zhang PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscTruth *flag)
1566*6e4ee0c6SHong Zhang {
1567*6e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
1568*6e4ee0c6SHong Zhang   Mat            a,b;
1569*6e4ee0c6SHong Zhang   PetscTruth     flg;
1570*6e4ee0c6SHong Zhang   PetscErrorCode ierr;
157190ace30eSBarry Smith 
1572*6e4ee0c6SHong Zhang   PetscFunctionBegin;
1573*6e4ee0c6SHong Zhang   a = matA->A;
1574*6e4ee0c6SHong Zhang   b = matB->A;
1575*6e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
1576*6e4ee0c6SHong Zhang   PetscFunctionReturn(0);
1577*6e4ee0c6SHong Zhang }
157890ace30eSBarry Smith 
157990ace30eSBarry Smith 
158090ace30eSBarry Smith 
1581