1*ca44d042SBarry Smith /*$Id: mpidense.c,v 1.139 2000/05/05 22:15:35 balay Exp bsmith $*/ 28965ea79SLois Curfman McInnes 3ed3cc1f0SBarry Smith /* 4ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 5ed3cc1f0SBarry Smith */ 6ed3cc1f0SBarry Smith 770f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 8f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 98965ea79SLois Curfman McInnes 100de54da6SSatish Balay EXTERN_C_BEGIN 110de54da6SSatish Balay #undef __FUNC__ 12b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetDiagonalBlock_MPIDense" 130de54da6SSatish Balay int MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 140de54da6SSatish Balay { 150de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 160de54da6SSatish Balay int m = mdn->m,rstart = mdn->rstart,rank,ierr; 170de54da6SSatish Balay Scalar *array; 180de54da6SSatish Balay MPI_Comm comm; 190de54da6SSatish Balay 200de54da6SSatish Balay PetscFunctionBegin; 210de54da6SSatish Balay if (mdn->M != mdn->N) SETERRQ(PETSC_ERR_SUP,0,"Only square matrices supported."); 220de54da6SSatish Balay 230de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 240de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 250de54da6SSatish Balay 260de54da6SSatish Balay ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 270de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 280de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 290de54da6SSatish Balay ierr = MatCreateSeqDense(comm,m,m,array+m*rstart,B);CHKERRQ(ierr); 300de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 310de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 320de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 330de54da6SSatish Balay 340de54da6SSatish Balay *iscopy = PETSC_TRUE; 350de54da6SSatish Balay PetscFunctionReturn(0); 360de54da6SSatish Balay } 370de54da6SSatish Balay EXTERN_C_END 380de54da6SSatish Balay 39*ca44d042SBarry Smith EXTERN int MatSetUpMultiply_MPIDense(Mat); 407ef1d9bdSSatish Balay 415615d1e5SSatish Balay #undef __FUNC__ 42b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValues_MPIDense" 438f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 448965ea79SLois Curfman McInnes { 4539b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 4639b7565bSBarry Smith int ierr,i,j,rstart = A->rstart,rend = A->rend,row; 4739b7565bSBarry Smith int roworiented = A->roworiented; 488965ea79SLois Curfman McInnes 493a40ed3dSBarry Smith PetscFunctionBegin; 508965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 515ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 52a8c6a408SBarry Smith if (idxm[i] >= A->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 538965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 548965ea79SLois Curfman McInnes row = idxm[i] - rstart; 5539b7565bSBarry Smith if (roworiented) { 5639b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 573a40ed3dSBarry Smith } else { 588965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 595ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 60a8c6a408SBarry Smith if (idxn[j] >= A->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6139b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 6239b7565bSBarry Smith } 638965ea79SLois Curfman McInnes } 643a40ed3dSBarry Smith } else { 653782ba37SSatish Balay if (!A->donotstash) { 6639b7565bSBarry Smith if (roworiented) { 678798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 68d36fbae8SSatish Balay } else { 698798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 7039b7565bSBarry Smith } 71b49de8d1SLois Curfman McInnes } 72b49de8d1SLois Curfman McInnes } 733782ba37SSatish Balay } 743a40ed3dSBarry Smith PetscFunctionReturn(0); 75b49de8d1SLois Curfman McInnes } 76b49de8d1SLois Curfman McInnes 775615d1e5SSatish Balay #undef __FUNC__ 78b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetValues_MPIDense" 798f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 80b49de8d1SLois Curfman McInnes { 81b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 82b49de8d1SLois Curfman McInnes int ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row; 83b49de8d1SLois Curfman McInnes 843a40ed3dSBarry Smith PetscFunctionBegin; 85b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 86a8c6a408SBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 87a8c6a408SBarry Smith if (idxm[i] >= mdn->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 88b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 89b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 90b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 91a8c6a408SBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 92a8c6a408SBarry Smith if (idxn[j] >= mdn->N) { 93a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 94a8c6a408SBarry Smith } 95b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 96b49de8d1SLois Curfman McInnes } 97a8c6a408SBarry Smith } else { 98a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported"); 998965ea79SLois Curfman McInnes } 1008965ea79SLois Curfman McInnes } 1013a40ed3dSBarry Smith PetscFunctionReturn(0); 1028965ea79SLois Curfman McInnes } 1038965ea79SLois Curfman McInnes 1045615d1e5SSatish Balay #undef __FUNC__ 105b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetArray_MPIDense" 1068f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array) 107ff14e315SSatish Balay { 108ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 109ff14e315SSatish Balay int ierr; 110ff14e315SSatish Balay 1113a40ed3dSBarry Smith PetscFunctionBegin; 112ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1133a40ed3dSBarry Smith PetscFunctionReturn(0); 114ff14e315SSatish Balay } 115ff14e315SSatish Balay 1165615d1e5SSatish Balay #undef __FUNC__ 117b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrix_MPIDense" 118ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 119ca3fa75bSLois Curfman McInnes { 120ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 121ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 12272d926a5SLois Curfman McInnes int i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols,rank; 123ca3fa75bSLois Curfman McInnes Scalar *av,*bv,*v = lmat->v; 124ca3fa75bSLois Curfman McInnes Mat newmat; 125ca3fa75bSLois Curfman McInnes 126ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 1277eba5e9cSLois Curfman McInnes ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 128ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 129ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 130ca3fa75bSLois Curfman McInnes ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr); 131ca3fa75bSLois Curfman McInnes ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr); 132ca3fa75bSLois Curfman McInnes 133ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1347eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1357eba5e9cSLois Curfman McInnes original matrix! */ 136ca3fa75bSLois Curfman McInnes 137ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1387eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 139ca3fa75bSLois Curfman McInnes 140ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 141ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 1427eba5e9cSLois Curfman McInnes /* SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); */ 1437eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 144ca3fa75bSLois Curfman McInnes newmat = *B; 145ca3fa75bSLois Curfman McInnes } else { 146ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 147ca3fa75bSLois Curfman McInnes ierr = MatCreateMPIDense(A->comm,nrows,ncols,PETSC_DECIDE,PETSC_DECIDE,PETSC_NULL,&newmat);CHKERRQ(ierr); 148ca3fa75bSLois Curfman McInnes } 149ca3fa75bSLois Curfman McInnes 150ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 151ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 152ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 153ca3fa75bSLois Curfman McInnes 154ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 155ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 156ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 1577eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 158ca3fa75bSLois Curfman McInnes } 159ca3fa75bSLois Curfman McInnes } 160ca3fa75bSLois Curfman McInnes 161ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 162ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164ca3fa75bSLois Curfman McInnes 165ca3fa75bSLois Curfman McInnes /* Free work space */ 166ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 167ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 168ca3fa75bSLois Curfman McInnes *B = newmat; 169ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 170ca3fa75bSLois Curfman McInnes } 171ca3fa75bSLois Curfman McInnes 172ca3fa75bSLois Curfman McInnes #undef __FUNC__ 173b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreArray_MPIDense" 1748f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array) 175ff14e315SSatish Balay { 1763a40ed3dSBarry Smith PetscFunctionBegin; 1773a40ed3dSBarry Smith PetscFunctionReturn(0); 178ff14e315SSatish Balay } 179ff14e315SSatish Balay 1805615d1e5SSatish Balay #undef __FUNC__ 181b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatAssemblyBegin_MPIDense" 1828f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 1838965ea79SLois Curfman McInnes { 18439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 1858965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 186d36fbae8SSatish Balay int ierr,nstash,reallocs; 1878965ea79SLois Curfman McInnes InsertMode addv; 1888965ea79SLois Curfman McInnes 1893a40ed3dSBarry Smith PetscFunctionBegin; 1908965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 191ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1927056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 193a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Cannot mix adds/inserts on different procs"); 1948965ea79SLois Curfman McInnes } 195e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1968965ea79SLois Curfman McInnes 1978798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 1988798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 199c38d4ed2SBarry Smith PLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs); 2003a40ed3dSBarry Smith PetscFunctionReturn(0); 2018965ea79SLois Curfman McInnes } 2028965ea79SLois Curfman McInnes 2035615d1e5SSatish Balay #undef __FUNC__ 204b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatAssemblyEnd_MPIDense" 2058f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2068965ea79SLois Curfman McInnes { 20739ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2087ef1d9bdSSatish Balay int i,n,ierr,*row,*col,flg,j,rstart,ncols; 2097ef1d9bdSSatish Balay Scalar *val; 210e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2118965ea79SLois Curfman McInnes 2123a40ed3dSBarry Smith PetscFunctionBegin; 2138965ea79SLois Curfman McInnes /* wait on receives */ 2147ef1d9bdSSatish Balay while (1) { 2158798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2167ef1d9bdSSatish Balay if (!flg) break; 2178965ea79SLois Curfman McInnes 2187ef1d9bdSSatish Balay for (i=0; i<n;) { 2197ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2207ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2217ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2227ef1d9bdSSatish Balay else ncols = n-i; 2237ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2247ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2257ef1d9bdSSatish Balay i = j; 2268965ea79SLois Curfman McInnes } 2277ef1d9bdSSatish Balay } 2288798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2298965ea79SLois Curfman McInnes 23039ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 23139ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2328965ea79SLois Curfman McInnes 2336d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 23439ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2358965ea79SLois Curfman McInnes } 2363a40ed3dSBarry Smith PetscFunctionReturn(0); 2378965ea79SLois Curfman McInnes } 2388965ea79SLois Curfman McInnes 2395615d1e5SSatish Balay #undef __FUNC__ 240b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroEntries_MPIDense" 2418f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A) 2428965ea79SLois Curfman McInnes { 2433a40ed3dSBarry Smith int ierr; 24439ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2453a40ed3dSBarry Smith 2463a40ed3dSBarry Smith PetscFunctionBegin; 2473a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2483a40ed3dSBarry Smith PetscFunctionReturn(0); 2498965ea79SLois Curfman McInnes } 2508965ea79SLois Curfman McInnes 2515615d1e5SSatish Balay #undef __FUNC__ 252b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetBlockSize_MPIDense" 2538f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs) 2544e220ebcSLois Curfman McInnes { 2553a40ed3dSBarry Smith PetscFunctionBegin; 2564e220ebcSLois Curfman McInnes *bs = 1; 2573a40ed3dSBarry Smith PetscFunctionReturn(0); 2584e220ebcSLois Curfman McInnes } 2594e220ebcSLois Curfman McInnes 2608965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2618965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2628965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2633501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2648965ea79SLois Curfman McInnes routine. 2658965ea79SLois Curfman McInnes */ 2665615d1e5SSatish Balay #undef __FUNC__ 267b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_MPIDense" 2688f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 2698965ea79SLois Curfman McInnes { 27039ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2718965ea79SLois Curfman McInnes int i,ierr,N,*rows,*owners = l->rowners,size = l->size; 2728965ea79SLois Curfman McInnes int *procs,*nprocs,j,found,idx,nsends,*work; 2738965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2748965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2758965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2768965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2778965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2788965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2798965ea79SLois Curfman McInnes IS istmp; 2808965ea79SLois Curfman McInnes 2813a40ed3dSBarry Smith PetscFunctionBegin; 28277c4ece6SBarry Smith ierr = ISGetSize(is,&N);CHKERRQ(ierr); 2838965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 2848965ea79SLois Curfman McInnes 2858965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 2860452661fSBarry Smith nprocs = (int*)PetscMalloc(2*size*sizeof(int));CHKPTRQ(nprocs); 287549d3d68SSatish Balay ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 288549d3d68SSatish Balay procs = nprocs + size; 2890452661fSBarry Smith owner = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(owner); /* see note*/ 2908965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 2918965ea79SLois Curfman McInnes idx = rows[i]; 2928965ea79SLois Curfman McInnes found = 0; 2938965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 2948965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 2958965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 2968965ea79SLois Curfman McInnes } 2978965ea79SLois Curfman McInnes } 298a8c6a408SBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range"); 2998965ea79SLois Curfman McInnes } 3008965ea79SLois Curfman McInnes nsends = 0; for (i=0; i<size; i++) { nsends += procs[i];} 3018965ea79SLois Curfman McInnes 3028965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 3036831982aSBarry Smith work = (int*)PetscMalloc(2*size*sizeof(int));CHKPTRQ(work); 3046831982aSBarry Smith ierr = MPI_Allreduce(nprocs,work,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr); 3058965ea79SLois Curfman McInnes nmax = work[rank]; 3066831982aSBarry Smith nrecvs = work[size+rank]; 307606d414cSSatish Balay ierr = PetscFree(work);CHKERRQ(ierr); 3088965ea79SLois Curfman McInnes 3098965ea79SLois Curfman McInnes /* post receives: */ 3103a40ed3dSBarry Smith rvalues = (int*)PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 3113a40ed3dSBarry Smith recv_waits = (MPI_Request*)PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 3128965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 313ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3148965ea79SLois Curfman McInnes } 3158965ea79SLois Curfman McInnes 3168965ea79SLois Curfman McInnes /* do sends: 3178965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3188965ea79SLois Curfman McInnes the ith processor 3198965ea79SLois Curfman McInnes */ 3200452661fSBarry Smith svalues = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(svalues); 3217056b6fcSBarry Smith send_waits = (MPI_Request*)PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 3220452661fSBarry Smith starts = (int*)PetscMalloc((size+1)*sizeof(int));CHKPTRQ(starts); 3238965ea79SLois Curfman McInnes starts[0] = 0; 3248965ea79SLois Curfman McInnes for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];} 3258965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3268965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3278965ea79SLois Curfman McInnes } 3288965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3298965ea79SLois Curfman McInnes 3308965ea79SLois Curfman McInnes starts[0] = 0; 3318965ea79SLois Curfman McInnes for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[i-1];} 3328965ea79SLois Curfman McInnes count = 0; 3338965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 3348965ea79SLois Curfman McInnes if (procs[i]) { 335ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3368965ea79SLois Curfman McInnes } 3378965ea79SLois Curfman McInnes } 338606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3398965ea79SLois Curfman McInnes 3408965ea79SLois Curfman McInnes base = owners[rank]; 3418965ea79SLois Curfman McInnes 3428965ea79SLois Curfman McInnes /* wait on receives */ 3430452661fSBarry Smith lens = (int*)PetscMalloc(2*(nrecvs+1)*sizeof(int));CHKPTRQ(lens); 3448965ea79SLois Curfman McInnes source = lens + nrecvs; 3458965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3468965ea79SLois Curfman McInnes while (count) { 347ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3488965ea79SLois Curfman McInnes /* unpack receives into our local space */ 349ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 3508965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3518965ea79SLois Curfman McInnes lens[imdex] = n; 3528965ea79SLois Curfman McInnes slen += n; 3538965ea79SLois Curfman McInnes count--; 3548965ea79SLois Curfman McInnes } 355606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3568965ea79SLois Curfman McInnes 3578965ea79SLois Curfman McInnes /* move the data into the send scatter */ 3580452661fSBarry Smith lrows = (int*)PetscMalloc((slen+1)*sizeof(int));CHKPTRQ(lrows); 3598965ea79SLois Curfman McInnes count = 0; 3608965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 3618965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3628965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 3638965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3648965ea79SLois Curfman McInnes } 3658965ea79SLois Curfman McInnes } 366606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 367606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 368606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 369606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3708965ea79SLois Curfman McInnes 3718965ea79SLois Curfman McInnes /* actually zap the local rows */ 372029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 3738965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 374606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3758965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 3768965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 3778965ea79SLois Curfman McInnes 3788965ea79SLois Curfman McInnes /* wait on sends */ 3798965ea79SLois Curfman McInnes if (nsends) { 3807056b6fcSBarry Smith send_status = (MPI_Status*)PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 381ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 382606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 3838965ea79SLois Curfman McInnes } 384606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 385606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 3868965ea79SLois Curfman McInnes 3873a40ed3dSBarry Smith PetscFunctionReturn(0); 3888965ea79SLois Curfman McInnes } 3898965ea79SLois Curfman McInnes 3905615d1e5SSatish Balay #undef __FUNC__ 391b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMult_MPIDense" 3928f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3938965ea79SLois Curfman McInnes { 39439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 3958965ea79SLois Curfman McInnes int ierr; 396c456f294SBarry Smith 3973a40ed3dSBarry Smith PetscFunctionBegin; 39843a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39943a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 40044cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4013a40ed3dSBarry Smith PetscFunctionReturn(0); 4028965ea79SLois Curfman McInnes } 4038965ea79SLois Curfman McInnes 4045615d1e5SSatish Balay #undef __FUNC__ 405b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultAdd_MPIDense" 4068f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4078965ea79SLois Curfman McInnes { 40839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4098965ea79SLois Curfman McInnes int ierr; 410c456f294SBarry Smith 4113a40ed3dSBarry Smith PetscFunctionBegin; 41243a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41343a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41444cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4153a40ed3dSBarry Smith PetscFunctionReturn(0); 4168965ea79SLois Curfman McInnes } 4178965ea79SLois Curfman McInnes 4185615d1e5SSatish Balay #undef __FUNC__ 419b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultTranspose_MPIDense" 4207c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 421096963f5SLois Curfman McInnes { 422096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 423096963f5SLois Curfman McInnes int ierr; 4243501a2bdSLois Curfman McInnes Scalar zero = 0.0; 425096963f5SLois Curfman McInnes 4263a40ed3dSBarry Smith PetscFunctionBegin; 4273501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 4287c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 429537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 430537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4313a40ed3dSBarry Smith PetscFunctionReturn(0); 432096963f5SLois Curfman McInnes } 433096963f5SLois Curfman McInnes 4345615d1e5SSatish Balay #undef __FUNC__ 435b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultTransposeAdd_MPIDense" 4367c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 437096963f5SLois Curfman McInnes { 438096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 439096963f5SLois Curfman McInnes int ierr; 440096963f5SLois Curfman McInnes 4413a40ed3dSBarry Smith PetscFunctionBegin; 4423501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4437c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 444537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 445537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4463a40ed3dSBarry Smith PetscFunctionReturn(0); 447096963f5SLois Curfman McInnes } 448096963f5SLois Curfman McInnes 4495615d1e5SSatish Balay #undef __FUNC__ 450b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetDiagonal_MPIDense" 4518f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 4528965ea79SLois Curfman McInnes { 45339ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 454096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 45544cd7ae7SLois Curfman McInnes int ierr,len,i,n,m = a->m,radd; 45644cd7ae7SLois Curfman McInnes Scalar *x,zero = 0.0; 457ed3cc1f0SBarry Smith 4583a40ed3dSBarry Smith PetscFunctionBegin; 45944cd7ae7SLois Curfman McInnes VecSet(&zero,v); 460096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x);CHKERRQ(ierr); 461096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 462a8c6a408SBarry Smith if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 46344cd7ae7SLois Curfman McInnes len = PetscMin(aloc->m,aloc->n); 4647ddc982cSLois Curfman McInnes radd = a->rstart*m; 46544cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 466096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 467096963f5SLois Curfman McInnes } 4689a8c540fSSatish Balay ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 4693a40ed3dSBarry Smith PetscFunctionReturn(0); 4708965ea79SLois Curfman McInnes } 4718965ea79SLois Curfman McInnes 4725615d1e5SSatish Balay #undef __FUNC__ 473b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDestroy_MPIDense" 474e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 4758965ea79SLois Curfman McInnes { 4763501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4778965ea79SLois Curfman McInnes int ierr; 478ed3cc1f0SBarry Smith 4793a40ed3dSBarry Smith PetscFunctionBegin; 48094d884c6SBarry Smith 48194d884c6SBarry Smith if (mat->mapping) { 48294d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 48394d884c6SBarry Smith } 48494d884c6SBarry Smith if (mat->bmapping) { 48594d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 48694d884c6SBarry Smith } 487aa482453SBarry Smith #if defined(PETSC_USE_LOG) 488e1311b90SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N); 4898965ea79SLois Curfman McInnes #endif 4908798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 491606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 4923501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 4933501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4943501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 495622d7880SLois Curfman McInnes if (mdn->factor) { 496606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 497606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 498606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 499606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 500622d7880SLois Curfman McInnes } 501606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 50261b13de0SBarry Smith if (mat->rmap) { 50361b13de0SBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 50461b13de0SBarry Smith } 50561b13de0SBarry Smith if (mat->cmap) { 50661b13de0SBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 50790f02eecSBarry Smith } 5088965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 5090452661fSBarry Smith PetscHeaderDestroy(mat); 5103a40ed3dSBarry Smith PetscFunctionReturn(0); 5118965ea79SLois Curfman McInnes } 51239ddd567SLois Curfman McInnes 5135615d1e5SSatish Balay #undef __FUNC__ 514b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense_Binary" 51539ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 5168965ea79SLois Curfman McInnes { 51739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 5188965ea79SLois Curfman McInnes int ierr; 5197056b6fcSBarry Smith 5203a40ed3dSBarry Smith PetscFunctionBegin; 52139ddd567SLois Curfman McInnes if (mdn->size == 1) { 52239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5238965ea79SLois Curfman McInnes } 524a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 5253a40ed3dSBarry Smith PetscFunctionReturn(0); 5268965ea79SLois Curfman McInnes } 5278965ea79SLois Curfman McInnes 5285615d1e5SSatish Balay #undef __FUNC__ 529b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense_ASCIIorDraworSocket" 530f1af5d2fSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,Viewer viewer) 5318965ea79SLois Curfman McInnes { 53239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 53377ed5343SBarry Smith int ierr,format,size = mdn->size,rank = mdn->rank; 53419bcc07fSBarry Smith ViewerType vtype; 535f1af5d2fSBarry Smith PetscTruth isascii,isdraw; 5368965ea79SLois Curfman McInnes 5373a40ed3dSBarry Smith PetscFunctionBegin; 538f1af5d2fSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 539f1af5d2fSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 540f1af5d2fSBarry Smith if (isascii) { 5413a40ed3dSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 542888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 543639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 5444e220ebcSLois Curfman McInnes MatInfo info; 545888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 5466831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m, 5476831982aSBarry Smith (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr); 5486831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 5493501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5503a40ed3dSBarry Smith PetscFunctionReturn(0); 55196f6c058SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO) { 5523a40ed3dSBarry Smith PetscFunctionReturn(0); 5538965ea79SLois Curfman McInnes } 554f1af5d2fSBarry Smith } else if (isdraw) { 555f1af5d2fSBarry Smith Draw draw; 556f1af5d2fSBarry Smith PetscTruth isnull; 557f1af5d2fSBarry Smith 558f1af5d2fSBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 559f1af5d2fSBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); 560f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 561f1af5d2fSBarry Smith } 56277ed5343SBarry Smith 5638965ea79SLois Curfman McInnes if (size == 1) { 56439ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5653a40ed3dSBarry Smith } else { 5668965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5678965ea79SLois Curfman McInnes Mat A; 56839ddd567SLois Curfman McInnes int M = mdn->M,N = mdn->N,m,row,i,nz,*cols; 56939ddd567SLois Curfman McInnes Scalar *vals; 57039ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*)mdn->A->data; 5718965ea79SLois Curfman McInnes 5728965ea79SLois Curfman McInnes if (!rank) { 573f1af5d2fSBarry Smith ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 5743a40ed3dSBarry Smith } else { 575f1af5d2fSBarry Smith ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 5768965ea79SLois Curfman McInnes } 5778965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 5788965ea79SLois Curfman McInnes 57939ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 58039ddd567SLois Curfman McInnes but it's quick for now */ 58139ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 5828965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 58339ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 58439ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 58539ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 58639ddd567SLois Curfman McInnes row++; 5878965ea79SLois Curfman McInnes } 5888965ea79SLois Curfman McInnes 5896d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5906d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5918965ea79SLois Curfman McInnes if (!rank) { 5926831982aSBarry Smith Viewer sviewer; 5936831982aSBarry Smith ierr = ViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 5946831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 5956831982aSBarry Smith ierr = ViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 5968965ea79SLois Curfman McInnes } 5976831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 5988965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 5998965ea79SLois Curfman McInnes } 6003a40ed3dSBarry Smith PetscFunctionReturn(0); 6018965ea79SLois Curfman McInnes } 6028965ea79SLois Curfman McInnes 6035615d1e5SSatish Balay #undef __FUNC__ 604b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIDense" 605e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer) 6068965ea79SLois Curfman McInnes { 60739ddd567SLois Curfman McInnes int ierr; 608f1af5d2fSBarry Smith PetscTruth isascii,isbinary,isdraw,issocket; 6098965ea79SLois Curfman McInnes 610433994e6SBarry Smith PetscFunctionBegin; 6110f5bd95cSBarry Smith 6126831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 6136831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 614f1af5d2fSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 615f1af5d2fSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 6160f5bd95cSBarry Smith 617f1af5d2fSBarry Smith if (isascii || issocket || isdraw) { 618f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 6190f5bd95cSBarry Smith } else if (isbinary) { 6203a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6215cd90555SBarry Smith } else { 6220f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6238965ea79SLois Curfman McInnes } 6243a40ed3dSBarry Smith PetscFunctionReturn(0); 6258965ea79SLois Curfman McInnes } 6268965ea79SLois Curfman McInnes 6275615d1e5SSatish Balay #undef __FUNC__ 628b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetInfo_MPIDense" 6298f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6308965ea79SLois Curfman McInnes { 6313501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 6323501a2bdSLois Curfman McInnes Mat mdn = mat->A; 6334e220ebcSLois Curfman McInnes int ierr; 634329f5518SBarry Smith PetscReal isend[5],irecv[5]; 6358965ea79SLois Curfman McInnes 6363a40ed3dSBarry Smith PetscFunctionBegin; 6374e220ebcSLois Curfman McInnes info->rows_global = (double)mat->M; 6384e220ebcSLois Curfman McInnes info->columns_global = (double)mat->N; 6394e220ebcSLois Curfman McInnes info->rows_local = (double)mat->m; 6404e220ebcSLois Curfman McInnes info->columns_local = (double)mat->N; 6414e220ebcSLois Curfman McInnes info->block_size = 1.0; 6424e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6434e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6444e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6458965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6464e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6474e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6484e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6494e220ebcSLois Curfman McInnes info->memory = isend[3]; 6504e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6518965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 652f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 6534e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6544e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6554e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6564e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6574e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6588965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 659f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 6604e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6614e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6624e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6634e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6644e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6658965ea79SLois Curfman McInnes } 6664e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6674e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6684e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6693a40ed3dSBarry Smith PetscFunctionReturn(0); 6708965ea79SLois Curfman McInnes } 6718965ea79SLois Curfman McInnes 6725615d1e5SSatish Balay #undef __FUNC__ 673b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_MPIDense" 6748f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 6758965ea79SLois Curfman McInnes { 67639ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 6778965ea79SLois Curfman McInnes 6783a40ed3dSBarry Smith PetscFunctionBegin; 6796d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 6806d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 6814787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 6824787f768SSatish Balay op == MAT_NEW_NONZERO_ALLOCATION_ERR || 683219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 684219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED) { 685b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 686b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 687aeafbbfcSLois Curfman McInnes a->roworiented = 1; 6888965ea79SLois Curfman McInnes MatSetOption(a->A,op); 689b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 690219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 6916d4a8577SBarry Smith op == MAT_SYMMETRIC || 6926d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 693b51ba29fSSatish Balay op == MAT_YES_NEW_DIAGONALS || 694b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 695981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 6963a40ed3dSBarry Smith } else if (op == MAT_COLUMN_ORIENTED) { 6973a40ed3dSBarry Smith a->roworiented = 0; MatSetOption(a->A,op); 6983782ba37SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 6993782ba37SSatish Balay a->donotstash = 1; 7003a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS) { 7013a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 7023a40ed3dSBarry Smith } else { 7033a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 7043a40ed3dSBarry Smith } 7053a40ed3dSBarry Smith PetscFunctionReturn(0); 7068965ea79SLois Curfman McInnes } 7078965ea79SLois Curfman McInnes 7085615d1e5SSatish Balay #undef __FUNC__ 709b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSize_MPIDense" 7108f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n) 7118965ea79SLois Curfman McInnes { 7123501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7133a40ed3dSBarry Smith 7143a40ed3dSBarry Smith PetscFunctionBegin; 715f1af5d2fSBarry Smith if (m) *m = mat->M; 716f1af5d2fSBarry Smith if (n) *n = mat->N; 7173a40ed3dSBarry Smith PetscFunctionReturn(0); 7188965ea79SLois Curfman McInnes } 7198965ea79SLois Curfman McInnes 7205615d1e5SSatish Balay #undef __FUNC__ 721b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetLocalSize_MPIDense" 7228f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 7238965ea79SLois Curfman McInnes { 7243501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7253a40ed3dSBarry Smith 7263a40ed3dSBarry Smith PetscFunctionBegin; 7278965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 7283a40ed3dSBarry Smith PetscFunctionReturn(0); 7298965ea79SLois Curfman McInnes } 7308965ea79SLois Curfman McInnes 7315615d1e5SSatish Balay #undef __FUNC__ 732b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_MPIDense" 7338f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 7348965ea79SLois Curfman McInnes { 7353501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7363a40ed3dSBarry Smith 7373a40ed3dSBarry Smith PetscFunctionBegin; 7384c49b128SBarry Smith if (m) *m = mat->rstart; 7394c49b128SBarry Smith if (n) *n = mat->rend; 7403a40ed3dSBarry Smith PetscFunctionReturn(0); 7418965ea79SLois Curfman McInnes } 7428965ea79SLois Curfman McInnes 7435615d1e5SSatish Balay #undef __FUNC__ 744b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_MPIDense" 7458f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 7468965ea79SLois Curfman McInnes { 7473501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7483a40ed3dSBarry Smith int lrow,rstart = mat->rstart,rend = mat->rend,ierr; 7498965ea79SLois Curfman McInnes 7503a40ed3dSBarry Smith PetscFunctionBegin; 751a8c6a408SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows") 7528965ea79SLois Curfman McInnes lrow = row - rstart; 7533a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 7543a40ed3dSBarry Smith PetscFunctionReturn(0); 7558965ea79SLois Curfman McInnes } 7568965ea79SLois Curfman McInnes 7575615d1e5SSatish Balay #undef __FUNC__ 758b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_MPIDense" 7598f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 7608965ea79SLois Curfman McInnes { 761606d414cSSatish Balay int ierr; 762606d414cSSatish Balay 7633a40ed3dSBarry Smith PetscFunctionBegin; 764606d414cSSatish Balay if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 765606d414cSSatish Balay if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 7663a40ed3dSBarry Smith PetscFunctionReturn(0); 7678965ea79SLois Curfman McInnes } 7688965ea79SLois Curfman McInnes 7695615d1e5SSatish Balay #undef __FUNC__ 770b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDiagonalScale_MPIDense" 7715b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7725b2fa520SLois Curfman McInnes { 7735b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7745b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 7755b2fa520SLois Curfman McInnes Scalar *l,*r,x,*v; 77672d926a5SLois Curfman McInnes int ierr,i,j,s2a,s3a,s2,s3,m=mat->m,n=mat->n; 7775b2fa520SLois Curfman McInnes 7785b2fa520SLois Curfman McInnes PetscFunctionBegin; 77972d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7805b2fa520SLois Curfman McInnes if (ll) { 78172d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 78272d926a5SLois Curfman McInnes if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector non-conforming local size"); 7835b2fa520SLois Curfman McInnes ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7845b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7855b2fa520SLois Curfman McInnes x = l[i]; 7865b2fa520SLois Curfman McInnes v = mat->v + i; 7875b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7885b2fa520SLois Curfman McInnes } 7895b2fa520SLois Curfman McInnes ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 7905b2fa520SLois Curfman McInnes PLogFlops(n*m); 7915b2fa520SLois Curfman McInnes } 7925b2fa520SLois Curfman McInnes if (rr) { 79372d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 79472d926a5SLois Curfman McInnes if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec non-conforming local size"); 7955b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7965b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7975b2fa520SLois Curfman McInnes ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 7985b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 7995b2fa520SLois Curfman McInnes x = r[i]; 8005b2fa520SLois Curfman McInnes v = mat->v + i*m; 8015b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 8025b2fa520SLois Curfman McInnes } 80372d926a5SLois Curfman McInnes ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 8045b2fa520SLois Curfman McInnes PLogFlops(n*m); 8055b2fa520SLois Curfman McInnes } 8065b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 8075b2fa520SLois Curfman McInnes } 8085b2fa520SLois Curfman McInnes 8095b2fa520SLois Curfman McInnes #undef __FUNC__ 810b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatNorm_MPIDense" 811329f5518SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *norm) 812096963f5SLois Curfman McInnes { 8133501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8143501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 8153501a2bdSLois Curfman McInnes int ierr,i,j; 816329f5518SBarry Smith PetscReal sum = 0.0; 8173501a2bdSLois Curfman McInnes Scalar *v = mat->v; 8183501a2bdSLois Curfman McInnes 8193a40ed3dSBarry Smith PetscFunctionBegin; 8203501a2bdSLois Curfman McInnes if (mdn->size == 1) { 8213501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm);CHKERRQ(ierr); 8223501a2bdSLois Curfman McInnes } else { 8233501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 8243501a2bdSLois Curfman McInnes for (i=0; i<mat->n*mat->m; i++) { 825aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 826329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 8273501a2bdSLois Curfman McInnes #else 8283501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 8293501a2bdSLois Curfman McInnes #endif 8303501a2bdSLois Curfman McInnes } 831ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 8323501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 8333501a2bdSLois Curfman McInnes PLogFlops(2*mat->n*mat->m); 8343a40ed3dSBarry Smith } else if (type == NORM_1) { 835329f5518SBarry Smith PetscReal *tmp,*tmp2; 836329f5518SBarry Smith tmp = (PetscReal*)PetscMalloc(2*mdn->N*sizeof(PetscReal));CHKPTRQ(tmp); 8373501a2bdSLois Curfman McInnes tmp2 = tmp + mdn->N; 838329f5518SBarry Smith ierr = PetscMemzero(tmp,2*mdn->N*sizeof(PetscReal));CHKERRQ(ierr); 839096963f5SLois Curfman McInnes *norm = 0.0; 8403501a2bdSLois Curfman McInnes v = mat->v; 8413501a2bdSLois Curfman McInnes for (j=0; j<mat->n; j++) { 8423501a2bdSLois Curfman McInnes for (i=0; i<mat->m; i++) { 84367e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8443501a2bdSLois Curfman McInnes } 8453501a2bdSLois Curfman McInnes } 846ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 8473501a2bdSLois Curfman McInnes for (j=0; j<mdn->N; j++) { 8483501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 8493501a2bdSLois Curfman McInnes } 850606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 8513501a2bdSLois Curfman McInnes PLogFlops(mat->n*mat->m); 8523a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 853329f5518SBarry Smith PetscReal ntemp; 8543501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 855ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 8563a40ed3dSBarry Smith } else { 857a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 8583501a2bdSLois Curfman McInnes } 8593501a2bdSLois Curfman McInnes } 8603a40ed3dSBarry Smith PetscFunctionReturn(0); 8613501a2bdSLois Curfman McInnes } 8623501a2bdSLois Curfman McInnes 8635615d1e5SSatish Balay #undef __FUNC__ 864b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatTranspose_MPIDense" 8658f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 8663501a2bdSLois Curfman McInnes { 8673501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8683501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8693501a2bdSLois Curfman McInnes Mat B; 8703501a2bdSLois Curfman McInnes int M = a->M,N = a->N,m,n,*rwork,rstart = a->rstart; 8713501a2bdSLois Curfman McInnes int j,i,ierr; 8723501a2bdSLois Curfman McInnes Scalar *v; 8733501a2bdSLois Curfman McInnes 8743a40ed3dSBarry Smith PetscFunctionBegin; 8757c922b88SBarry Smith if (!matout && M != N) { 876a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place"); 8777056b6fcSBarry Smith } 8787056b6fcSBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 8793501a2bdSLois Curfman McInnes 8803501a2bdSLois Curfman McInnes m = Aloc->m; n = Aloc->n; v = Aloc->v; 8810452661fSBarry Smith rwork = (int*)PetscMalloc(n*sizeof(int));CHKPTRQ(rwork); 8823501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 8833501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8843501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8853501a2bdSLois Curfman McInnes v += m; 8863501a2bdSLois Curfman McInnes } 887606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 8886d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8896d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8907c922b88SBarry Smith if (matout) { 8913501a2bdSLois Curfman McInnes *matout = B; 8923501a2bdSLois Curfman McInnes } else { 893f830108cSBarry Smith PetscOps *Abops; 89409dc0095SBarry Smith MatOps Aops; 895f830108cSBarry Smith 8963501a2bdSLois Curfman McInnes /* This isn't really an in-place transpose, but free data struct from a */ 897606d414cSSatish Balay ierr = PetscFree(a->rowners);CHKERRQ(ierr); 8983501a2bdSLois Curfman McInnes ierr = MatDestroy(a->A);CHKERRQ(ierr); 8993501a2bdSLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 9003501a2bdSLois Curfman McInnes if (a->Mvctx) VecScatterDestroy(a->Mvctx); 901606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 902f830108cSBarry Smith 903f830108cSBarry Smith /* 904f830108cSBarry Smith This is horrible, horrible code. We need to keep the 905f830108cSBarry Smith A pointers for the bops and ops but copy everything 906f830108cSBarry Smith else from C. 907f830108cSBarry Smith */ 908f830108cSBarry Smith Abops = A->bops; 909f830108cSBarry Smith Aops = A->ops; 910549d3d68SSatish Balay ierr = PetscMemcpy(A,B,sizeof(struct _p_Mat));CHKERRQ(ierr); 911f830108cSBarry Smith A->bops = Abops; 912f830108cSBarry Smith A->ops = Aops; 913f830108cSBarry Smith 9140452661fSBarry Smith PetscHeaderDestroy(B); 9153501a2bdSLois Curfman McInnes } 9163a40ed3dSBarry Smith PetscFunctionReturn(0); 917096963f5SLois Curfman McInnes } 918096963f5SLois Curfman McInnes 919eadb2fb4SBarry Smith #include "pinclude/blaslapack.h" 9205615d1e5SSatish Balay #undef __FUNC__ 921b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatScale_MPIDense" 9228f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA) 92344cd7ae7SLois Curfman McInnes { 92444cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 92544cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 92644cd7ae7SLois Curfman McInnes int one = 1,nz; 92744cd7ae7SLois Curfman McInnes 9283a40ed3dSBarry Smith PetscFunctionBegin; 92944cd7ae7SLois Curfman McInnes nz = a->m*a->n; 93044cd7ae7SLois Curfman McInnes BLscal_(&nz,alpha,a->v,&one); 93144cd7ae7SLois Curfman McInnes PLogFlops(nz); 9323a40ed3dSBarry Smith PetscFunctionReturn(0); 93344cd7ae7SLois Curfman McInnes } 93444cd7ae7SLois Curfman McInnes 9355609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 936*ca44d042SBarry Smith EXTERN int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **); 9378965ea79SLois Curfman McInnes 9388965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 93909dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 94009dc0095SBarry Smith MatGetRow_MPIDense, 94109dc0095SBarry Smith MatRestoreRow_MPIDense, 94209dc0095SBarry Smith MatMult_MPIDense, 94309dc0095SBarry Smith MatMultAdd_MPIDense, 9447c922b88SBarry Smith MatMultTranspose_MPIDense, 9457c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 9468965ea79SLois Curfman McInnes 0, 94709dc0095SBarry Smith 0, 94809dc0095SBarry Smith 0, 94909dc0095SBarry Smith 0, 95009dc0095SBarry Smith 0, 95109dc0095SBarry Smith 0, 95209dc0095SBarry Smith 0, 95309dc0095SBarry Smith MatTranspose_MPIDense, 95409dc0095SBarry Smith MatGetInfo_MPIDense,0, 95509dc0095SBarry Smith MatGetDiagonal_MPIDense, 9565b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 95709dc0095SBarry Smith MatNorm_MPIDense, 95809dc0095SBarry Smith MatAssemblyBegin_MPIDense, 95909dc0095SBarry Smith MatAssemblyEnd_MPIDense, 96009dc0095SBarry Smith 0, 96109dc0095SBarry Smith MatSetOption_MPIDense, 96209dc0095SBarry Smith MatZeroEntries_MPIDense, 96309dc0095SBarry Smith MatZeroRows_MPIDense, 96409dc0095SBarry Smith 0, 96509dc0095SBarry Smith 0, 96609dc0095SBarry Smith 0, 96709dc0095SBarry Smith 0, 96809dc0095SBarry Smith MatGetSize_MPIDense, 96909dc0095SBarry Smith MatGetLocalSize_MPIDense, 97039ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 97109dc0095SBarry Smith 0, 97209dc0095SBarry Smith 0, 97309dc0095SBarry Smith MatGetArray_MPIDense, 97409dc0095SBarry Smith MatRestoreArray_MPIDense, 9755609ef8eSBarry Smith MatDuplicate_MPIDense, 97609dc0095SBarry Smith 0, 97709dc0095SBarry Smith 0, 97809dc0095SBarry Smith 0, 97909dc0095SBarry Smith 0, 98009dc0095SBarry Smith 0, 9812ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 98209dc0095SBarry Smith 0, 98309dc0095SBarry Smith MatGetValues_MPIDense, 98409dc0095SBarry Smith 0, 98509dc0095SBarry Smith 0, 98609dc0095SBarry Smith MatScale_MPIDense, 98709dc0095SBarry Smith 0, 98809dc0095SBarry Smith 0, 98909dc0095SBarry Smith 0, 99009dc0095SBarry Smith MatGetBlockSize_MPIDense, 99109dc0095SBarry Smith 0, 99209dc0095SBarry Smith 0, 99309dc0095SBarry Smith 0, 99409dc0095SBarry Smith 0, 99509dc0095SBarry Smith 0, 99609dc0095SBarry Smith 0, 99709dc0095SBarry Smith 0, 99809dc0095SBarry Smith 0, 99909dc0095SBarry Smith 0, 1000ca3fa75bSLois Curfman McInnes MatGetSubMatrix_MPIDense, 100109dc0095SBarry Smith 0, 100209dc0095SBarry Smith 0, 100309dc0095SBarry Smith MatGetMaps_Petsc}; 10048965ea79SLois Curfman McInnes 10055615d1e5SSatish Balay #undef __FUNC__ 1006b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateMPIDense" 10078965ea79SLois Curfman McInnes /*@C 100839ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 10098965ea79SLois Curfman McInnes 1010db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1011db81eaa0SLois Curfman McInnes 10128965ea79SLois Curfman McInnes Input Parameters: 1013db81eaa0SLois Curfman McInnes + comm - MPI communicator 10148965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1015db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 10168965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1017db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 1018db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1019dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 10208965ea79SLois Curfman McInnes 10218965ea79SLois Curfman McInnes Output Parameter: 1022477f1c0bSLois Curfman McInnes . A - the matrix 10238965ea79SLois Curfman McInnes 1024b259b22eSLois Curfman McInnes Notes: 102539ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 102639ddd567SLois Curfman McInnes storage by columns. 10278965ea79SLois Curfman McInnes 102818f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 102918f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1030b4fd4287SBarry Smith set data=PETSC_NULL. 103118f449edSLois Curfman McInnes 10328965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 10338965ea79SLois Curfman McInnes (possibly both). 10348965ea79SLois Curfman McInnes 1035027ccd11SLois Curfman McInnes Level: intermediate 1036027ccd11SLois Curfman McInnes 103739ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 10388965ea79SLois Curfman McInnes 103939ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 10408965ea79SLois Curfman McInnes @*/ 1041477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A) 10428965ea79SLois Curfman McInnes { 10438965ea79SLois Curfman McInnes Mat mat; 104439ddd567SLois Curfman McInnes Mat_MPIDense *a; 1045f1af5d2fSBarry Smith int ierr,i; 1046f1af5d2fSBarry Smith PetscTruth flg; 10478965ea79SLois Curfman McInnes 10483a40ed3dSBarry Smith PetscFunctionBegin; 1049ed2daf61SLois Curfman McInnes /* Note: For now, when data is specified above, this assumes the user correctly 1050ed2daf61SLois Curfman McInnes allocates the local dense storage space. We should add error checking. */ 105118f449edSLois Curfman McInnes 1052477f1c0bSLois Curfman McInnes *A = 0; 10533f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",comm,MatDestroy,MatView); 10548965ea79SLois Curfman McInnes PLogObjectCreate(mat); 10550452661fSBarry Smith mat->data = (void*)(a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 1056549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1057e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 1058e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 10598965ea79SLois Curfman McInnes mat->factor = 0; 106090f02eecSBarry Smith mat->mapping = 0; 10618965ea79SLois Curfman McInnes 1062622d7880SLois Curfman McInnes a->factor = 0; 1063e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1064d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&a->rank);CHKERRQ(ierr); 1065d132466eSBarry Smith ierr = MPI_Comm_size(comm,&a->size);CHKERRQ(ierr); 10668965ea79SLois Curfman McInnes 106796f6c058SBarry Smith ierr = PetscSplitOwnership(comm,&m,&M);CHKERRQ(ierr); 106839ddd567SLois Curfman McInnes 1069f1af5d2fSBarry Smith ierr = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr); 1070be0abb6dSBarry Smith a->nvec = n; 1071c7fcc2eaSBarry Smith 107239ddd567SLois Curfman McInnes /* each row stores all columns */ 1073aca0ad90SLois Curfman McInnes a->N = mat->N = N; 1074aca0ad90SLois Curfman McInnes a->M = mat->M = M; 1075aca0ad90SLois Curfman McInnes a->m = mat->m = m; 1076be0abb6dSBarry Smith a->n = mat->n = N; /* NOTE: n == N */ 10778965ea79SLois Curfman McInnes 1078c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 1079c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 1080488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,M,&mat->rmap);CHKERRQ(ierr); 1081be0abb6dSBarry Smith ierr = MapCreateMPI(comm,n,N,&mat->cmap);CHKERRQ(ierr); 1082c7fcc2eaSBarry Smith 10838965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 1084d7e8b826SBarry Smith a->rowners = (int*)PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners); 1085d7e8b826SBarry Smith a->cowners = a->rowners + a->size + 1; 1086f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1087ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 10888965ea79SLois Curfman McInnes a->rowners[0] = 0; 10898965ea79SLois Curfman McInnes for (i=2; i<=a->size; i++) { 10908965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 10918965ea79SLois Curfman McInnes } 10928965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 10938965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 1094ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1095d7e8b826SBarry Smith a->cowners[0] = 0; 1096d7e8b826SBarry Smith for (i=2; i<=a->size; i++) { 1097d7e8b826SBarry Smith a->cowners[i] += a->cowners[i-1]; 1098d7e8b826SBarry Smith } 10998965ea79SLois Curfman McInnes 1100029af93fSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A);CHKERRQ(ierr); 11018965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 11028965ea79SLois Curfman McInnes 11038965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 11043782ba37SSatish Balay a->donotstash = 0; 11058798bf22SSatish Balay ierr = MatStashCreate_Private(comm,1,&mat->stash);CHKERRQ(ierr); 11068965ea79SLois Curfman McInnes 11078965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 11088965ea79SLois Curfman McInnes a->lvec = 0; 11098965ea79SLois Curfman McInnes a->Mvctx = 0; 111039b7565bSBarry Smith a->roworiented = 1; 11118965ea79SLois Curfman McInnes 1112f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 11130de54da6SSatish Balay "MatGetDiagonalBlock_MPIDense", 11140c97f09dSSatish Balay MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 11150de54da6SSatish Balay 1116477f1c0bSLois Curfman McInnes *A = mat; 111725cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 111825cdf11fSBarry Smith if (flg) { 11198c469469SLois Curfman McInnes ierr = MatPrintHelp(mat);CHKERRQ(ierr); 11208c469469SLois Curfman McInnes } 11213a40ed3dSBarry Smith PetscFunctionReturn(0); 11228965ea79SLois Curfman McInnes } 11238965ea79SLois Curfman McInnes 11245615d1e5SSatish Balay #undef __FUNC__ 1125b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDuplicate_MPIDense" 11265609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 11278965ea79SLois Curfman McInnes { 11288965ea79SLois Curfman McInnes Mat mat; 11293501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 113039ddd567SLois Curfman McInnes int ierr; 11312ba99913SLois Curfman McInnes FactorCtx *factor; 11328965ea79SLois Curfman McInnes 11333a40ed3dSBarry Smith PetscFunctionBegin; 11348965ea79SLois Curfman McInnes *newmat = 0; 11353f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",A->comm,MatDestroy,MatView); 11368965ea79SLois Curfman McInnes PLogObjectCreate(mat); 11370452661fSBarry Smith mat->data = (void*)(a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 1138549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1139e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 1140e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 11413501a2bdSLois Curfman McInnes mat->factor = A->factor; 1142c456f294SBarry Smith mat->assembled = PETSC_TRUE; 11438965ea79SLois Curfman McInnes 114444cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 114544cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 114644cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 114744cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 11482ba99913SLois Curfman McInnes if (oldmat->factor) { 11492ba99913SLois Curfman McInnes a->factor = (FactorCtx*)(factor = PetscNew(FactorCtx));CHKPTRQ(factor); 11502ba99913SLois Curfman McInnes /* copy factor contents ... add this code! */ 11512ba99913SLois Curfman McInnes } else a->factor = 0; 11528965ea79SLois Curfman McInnes 11538965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 11548965ea79SLois Curfman McInnes a->rend = oldmat->rend; 11558965ea79SLois Curfman McInnes a->size = oldmat->size; 11568965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1157e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 11583782ba37SSatish Balay a->donotstash = oldmat->donotstash; 11590452661fSBarry Smith a->rowners = (int*)PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners); 1160f09e8eb9SSatish Balay PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1161549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 11628798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 11638965ea79SLois Curfman McInnes 1164329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 11655609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 11668965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 11678965ea79SLois Curfman McInnes *newmat = mat; 11683a40ed3dSBarry Smith PetscFunctionReturn(0); 11698965ea79SLois Curfman McInnes } 11708965ea79SLois Curfman McInnes 1171e090d566SSatish Balay #include "petscsys.h" 11728965ea79SLois Curfman McInnes 11735615d1e5SSatish Balay #undef __FUNC__ 1174b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_MPIDense_DenseInFile" 117590ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat) 117690ace30eSBarry Smith { 117740011551SBarry Smith int *rowners,i,size,rank,m,ierr,nz,j; 117890ace30eSBarry Smith Scalar *array,*vals,*vals_ptr; 117990ace30eSBarry Smith MPI_Status status; 118090ace30eSBarry Smith 11813a40ed3dSBarry Smith PetscFunctionBegin; 1182d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1183d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 118490ace30eSBarry Smith 118590ace30eSBarry Smith /* determine ownership of all rows */ 118690ace30eSBarry Smith m = M/size + ((M % size) > rank); 118790ace30eSBarry Smith rowners = (int*)PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1188ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 118990ace30eSBarry Smith rowners[0] = 0; 119090ace30eSBarry Smith for (i=2; i<=size; i++) { 119190ace30eSBarry Smith rowners[i] += rowners[i-1]; 119290ace30eSBarry Smith } 119390ace30eSBarry Smith 119490ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 119590ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 119690ace30eSBarry Smith 119790ace30eSBarry Smith if (!rank) { 119890ace30eSBarry Smith vals = (Scalar*)PetscMalloc(m*N*sizeof(Scalar));CHKPTRQ(vals); 119990ace30eSBarry Smith 120090ace30eSBarry Smith /* read in my part of the matrix numerical values */ 12010752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 120290ace30eSBarry Smith 120390ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 120490ace30eSBarry Smith vals_ptr = vals; 120590ace30eSBarry Smith for (i=0; i<m; i++) { 120690ace30eSBarry Smith for (j=0; j<N; j++) { 120790ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 120890ace30eSBarry Smith } 120990ace30eSBarry Smith } 121090ace30eSBarry Smith 121190ace30eSBarry Smith /* read in other processors and ship out */ 121290ace30eSBarry Smith for (i=1; i<size; i++) { 121390ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 12140752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1215ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 121690ace30eSBarry Smith } 12173a40ed3dSBarry Smith } else { 121890ace30eSBarry Smith /* receive numeric values */ 121990ace30eSBarry Smith vals = (Scalar*)PetscMalloc(m*N*sizeof(Scalar));CHKPTRQ(vals); 122090ace30eSBarry Smith 122190ace30eSBarry Smith /* receive message of values*/ 1222ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 122390ace30eSBarry Smith 122490ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 122590ace30eSBarry Smith vals_ptr = vals; 122690ace30eSBarry Smith for (i=0; i<m; i++) { 122790ace30eSBarry Smith for (j=0; j<N; j++) { 122890ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 122990ace30eSBarry Smith } 123090ace30eSBarry Smith } 123190ace30eSBarry Smith } 1232606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1233606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 12346d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12356d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12363a40ed3dSBarry Smith PetscFunctionReturn(0); 123790ace30eSBarry Smith } 123890ace30eSBarry Smith 123990ace30eSBarry Smith 12405615d1e5SSatish Balay #undef __FUNC__ 1241b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_MPIDense" 124219bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat) 12438965ea79SLois Curfman McInnes { 12448965ea79SLois Curfman McInnes Mat A; 12458965ea79SLois Curfman McInnes Scalar *vals,*svals; 124619bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 12478965ea79SLois Curfman McInnes MPI_Status status; 12488965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 12498965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols; 125019bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 12513a40ed3dSBarry Smith int i,nz,ierr,j,rstart,rend,fd; 12528965ea79SLois Curfman McInnes 12533a40ed3dSBarry Smith PetscFunctionBegin; 1254d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1255d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 12568965ea79SLois Curfman McInnes if (!rank) { 125790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 12580752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1259a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 12608965ea79SLois Curfman McInnes } 12618965ea79SLois Curfman McInnes 1262ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 126390ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 126490ace30eSBarry Smith 126590ace30eSBarry Smith /* 126690ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 126790ace30eSBarry Smith */ 126890ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 12693a40ed3dSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 12703a40ed3dSBarry Smith PetscFunctionReturn(0); 127190ace30eSBarry Smith } 127290ace30eSBarry Smith 12738965ea79SLois Curfman McInnes /* determine ownership of all rows */ 12748965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 12750452661fSBarry Smith rowners = (int*)PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1276ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 12778965ea79SLois Curfman McInnes rowners[0] = 0; 12788965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 12798965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 12808965ea79SLois Curfman McInnes } 12818965ea79SLois Curfman McInnes rstart = rowners[rank]; 12828965ea79SLois Curfman McInnes rend = rowners[rank+1]; 12838965ea79SLois Curfman McInnes 12848965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 12850452661fSBarry Smith ourlens = (int*)PetscMalloc(2*(rend-rstart)*sizeof(int));CHKPTRQ(ourlens); 12868965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 12878965ea79SLois Curfman McInnes if (!rank) { 12880452661fSBarry Smith rowlengths = (int*)PetscMalloc(M*sizeof(int));CHKPTRQ(rowlengths); 12890752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 12900452661fSBarry Smith sndcounts = (int*)PetscMalloc(size*sizeof(int));CHKPTRQ(sndcounts); 12918965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 1292ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1293606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1294ca161407SBarry Smith } else { 1295ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 12968965ea79SLois Curfman McInnes } 12978965ea79SLois Curfman McInnes 12988965ea79SLois Curfman McInnes if (!rank) { 12998965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 13000452661fSBarry Smith procsnz = (int*)PetscMalloc(size*sizeof(int));CHKPTRQ(procsnz); 1301549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 13028965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 13038965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 13048965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 13058965ea79SLois Curfman McInnes } 13068965ea79SLois Curfman McInnes } 1307606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 13088965ea79SLois Curfman McInnes 13098965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 13108965ea79SLois Curfman McInnes maxnz = 0; 13118965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 13120452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 13138965ea79SLois Curfman McInnes } 13140452661fSBarry Smith cols = (int*)PetscMalloc(maxnz*sizeof(int));CHKPTRQ(cols); 13158965ea79SLois Curfman McInnes 13168965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 13178965ea79SLois Curfman McInnes nz = procsnz[0]; 13180452661fSBarry Smith mycols = (int*)PetscMalloc(nz*sizeof(int));CHKPTRQ(mycols); 13190752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 13208965ea79SLois Curfman McInnes 13218965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 13228965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 13238965ea79SLois Curfman McInnes nz = procsnz[i]; 13240752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1325ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 13268965ea79SLois Curfman McInnes } 1327606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 13283a40ed3dSBarry Smith } else { 13298965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 13308965ea79SLois Curfman McInnes nz = 0; 13318965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13328965ea79SLois Curfman McInnes nz += ourlens[i]; 13338965ea79SLois Curfman McInnes } 13340452661fSBarry Smith mycols = (int*)PetscMalloc(nz*sizeof(int));CHKPTRQ(mycols); 13358965ea79SLois Curfman McInnes 13368965ea79SLois Curfman McInnes /* receive message of column indices*/ 1337ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1338ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1339a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 13408965ea79SLois Curfman McInnes } 13418965ea79SLois Curfman McInnes 13428965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1343549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 13448965ea79SLois Curfman McInnes jj = 0; 13458965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13468965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 13478965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 13488965ea79SLois Curfman McInnes jj++; 13498965ea79SLois Curfman McInnes } 13508965ea79SLois Curfman McInnes } 13518965ea79SLois Curfman McInnes 13528965ea79SLois Curfman McInnes /* create our matrix */ 13538965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13548965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 13558965ea79SLois Curfman McInnes } 1356b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 13578965ea79SLois Curfman McInnes A = *newmat; 13588965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13598965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 13608965ea79SLois Curfman McInnes } 13618965ea79SLois Curfman McInnes 13628965ea79SLois Curfman McInnes if (!rank) { 13630452661fSBarry Smith vals = (Scalar*)PetscMalloc(maxnz*sizeof(Scalar));CHKPTRQ(vals); 13648965ea79SLois Curfman McInnes 13658965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 13668965ea79SLois Curfman McInnes nz = procsnz[0]; 13670752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 13688965ea79SLois Curfman McInnes 13698965ea79SLois Curfman McInnes /* insert into matrix */ 13708965ea79SLois Curfman McInnes jj = rstart; 13718965ea79SLois Curfman McInnes smycols = mycols; 13728965ea79SLois Curfman McInnes svals = vals; 13738965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13748965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 13758965ea79SLois Curfman McInnes smycols += ourlens[i]; 13768965ea79SLois Curfman McInnes svals += ourlens[i]; 13778965ea79SLois Curfman McInnes jj++; 13788965ea79SLois Curfman McInnes } 13798965ea79SLois Curfman McInnes 13808965ea79SLois Curfman McInnes /* read in other processors and ship out */ 13818965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 13828965ea79SLois Curfman McInnes nz = procsnz[i]; 13830752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1384ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 13858965ea79SLois Curfman McInnes } 1386606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 13873a40ed3dSBarry Smith } else { 13888965ea79SLois Curfman McInnes /* receive numeric values */ 13890452661fSBarry Smith vals = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(vals); 13908965ea79SLois Curfman McInnes 13918965ea79SLois Curfman McInnes /* receive message of values*/ 1392ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1393ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1394a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 13958965ea79SLois Curfman McInnes 13968965ea79SLois Curfman McInnes /* insert into matrix */ 13978965ea79SLois Curfman McInnes jj = rstart; 13988965ea79SLois Curfman McInnes smycols = mycols; 13998965ea79SLois Curfman McInnes svals = vals; 14008965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14018965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14028965ea79SLois Curfman McInnes smycols += ourlens[i]; 14038965ea79SLois Curfman McInnes svals += ourlens[i]; 14048965ea79SLois Curfman McInnes jj++; 14058965ea79SLois Curfman McInnes } 14068965ea79SLois Curfman McInnes } 1407606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1408606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1409606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1410606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 14118965ea79SLois Curfman McInnes 14126d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14136d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14143a40ed3dSBarry Smith PetscFunctionReturn(0); 14158965ea79SLois Curfman McInnes } 141690ace30eSBarry Smith 141790ace30eSBarry Smith 141890ace30eSBarry Smith 141990ace30eSBarry Smith 142090ace30eSBarry Smith 1421