1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*7eba5e9cSLois Curfman McInnes static char vcid[] = "$Id: mpidense.c,v 1.119 1999/07/26 02:52:49 curfman Exp curfman $"; 38965ea79SLois Curfman McInnes #endif 48965ea79SLois Curfman McInnes 5ed3cc1f0SBarry Smith /* 6ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 7ed3cc1f0SBarry Smith */ 8ed3cc1f0SBarry Smith 970f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 10f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 118965ea79SLois Curfman McInnes 127ef1d9bdSSatish Balay extern int MatSetUpMultiply_MPIDense(Mat); 137ef1d9bdSSatish Balay 145615d1e5SSatish Balay #undef __FUNC__ 155615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIDense" 168f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 178965ea79SLois Curfman McInnes { 1839b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense *) mat->data; 1939b7565bSBarry Smith int ierr, i, j, rstart = A->rstart, rend = A->rend, row; 2039b7565bSBarry Smith int roworiented = A->roworiented; 218965ea79SLois Curfman McInnes 223a40ed3dSBarry Smith PetscFunctionBegin; 238965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 245ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 25a8c6a408SBarry Smith if (idxm[i] >= A->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 268965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 278965ea79SLois Curfman McInnes row = idxm[i] - rstart; 2839b7565bSBarry Smith if (roworiented) { 2939b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 303a40ed3dSBarry Smith } else { 318965ea79SLois Curfman McInnes for ( j=0; j<n; j++ ) { 325ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 33a8c6a408SBarry Smith if (idxn[j] >= A->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 3439b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 3539b7565bSBarry Smith } 368965ea79SLois Curfman McInnes } 373a40ed3dSBarry Smith } else { 383782ba37SSatish Balay if ( !A->donotstash) { 3939b7565bSBarry Smith if (roworiented) { 408798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 41d36fbae8SSatish Balay } else { 428798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 4339b7565bSBarry Smith } 44b49de8d1SLois Curfman McInnes } 45b49de8d1SLois Curfman McInnes } 463782ba37SSatish Balay } 473a40ed3dSBarry Smith PetscFunctionReturn(0); 48b49de8d1SLois Curfman McInnes } 49b49de8d1SLois Curfman McInnes 505615d1e5SSatish Balay #undef __FUNC__ 515615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIDense" 528f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 53b49de8d1SLois Curfman McInnes { 54b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 55b49de8d1SLois Curfman McInnes int ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row; 56b49de8d1SLois Curfman McInnes 573a40ed3dSBarry Smith PetscFunctionBegin; 58b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 59a8c6a408SBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 60a8c6a408SBarry Smith if (idxm[i] >= mdn->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 61b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 62b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 63b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 64a8c6a408SBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 65a8c6a408SBarry Smith if (idxn[j] >= mdn->N) { 66a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 67a8c6a408SBarry Smith } 68b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 69b49de8d1SLois Curfman McInnes } 70a8c6a408SBarry Smith } else { 71a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported"); 728965ea79SLois Curfman McInnes } 738965ea79SLois Curfman McInnes } 743a40ed3dSBarry Smith PetscFunctionReturn(0); 758965ea79SLois Curfman McInnes } 768965ea79SLois Curfman McInnes 775615d1e5SSatish Balay #undef __FUNC__ 785615d1e5SSatish Balay #define __FUNC__ "MatGetArray_MPIDense" 798f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array) 80ff14e315SSatish Balay { 81ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense *) A->data; 82ff14e315SSatish Balay int ierr; 83ff14e315SSatish Balay 843a40ed3dSBarry Smith PetscFunctionBegin; 85ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 863a40ed3dSBarry Smith PetscFunctionReturn(0); 87ff14e315SSatish Balay } 88ff14e315SSatish Balay 895615d1e5SSatish Balay #undef __FUNC__ 90ca3fa75bSLois Curfman McInnes #define __FUNC__ "MatGetSubMatrix_MPIDense" 91ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 92ca3fa75bSLois Curfman McInnes { 93ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data, *newmatd; 94ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense *) mat->A->data; 95*7eba5e9cSLois Curfman McInnes int i, j, ierr, *irow, *icol, n_lrows, n_lcols; 96*7eba5e9cSLois Curfman McInnes int rstart, rend, nrows, ncols, nlrows, nlcols, rank, ncols_all; 97ca3fa75bSLois Curfman McInnes Scalar *av, *bv, *v = lmat->v; 98ca3fa75bSLois Curfman McInnes Mat newmat; 99ca3fa75bSLois Curfman McInnes 100ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 101*7eba5e9cSLois Curfman McInnes ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 102ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 103ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 104ca3fa75bSLois Curfman McInnes ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr); 105ca3fa75bSLois Curfman McInnes ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr); 106ca3fa75bSLois Curfman McInnes 107ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 108*7eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 109*7eba5e9cSLois Curfman McInnes original matrix! */ 110ca3fa75bSLois Curfman McInnes 111ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 112*7eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 113ca3fa75bSLois Curfman McInnes 114ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 115ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 116*7eba5e9cSLois Curfman McInnes /* SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); */ 117*7eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 118ca3fa75bSLois Curfman McInnes newmat = *B; 119ca3fa75bSLois Curfman McInnes } else { 120ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 121ca3fa75bSLois Curfman McInnes ierr = MatCreateMPIDense(A->comm,nrows,ncols,PETSC_DECIDE,PETSC_DECIDE,PETSC_NULL,&newmat);CHKERRQ(ierr); 122ca3fa75bSLois Curfman McInnes } 123ca3fa75bSLois Curfman McInnes 124ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 125ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense *) newmat->data; 126ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 127ca3fa75bSLois Curfman McInnes 128ca3fa75bSLois Curfman McInnes for ( i=0; i<ncols; i++ ) { 129ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 130ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++ ) { 131*7eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 132ca3fa75bSLois Curfman McInnes } 133ca3fa75bSLois Curfman McInnes } 134ca3fa75bSLois Curfman McInnes 135ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 136ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 137ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 138ca3fa75bSLois Curfman McInnes 139ca3fa75bSLois Curfman McInnes /* Free work space */ 140ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 141ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 142ca3fa75bSLois Curfman McInnes *B = newmat; 143ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 144ca3fa75bSLois Curfman McInnes } 145ca3fa75bSLois Curfman McInnes 146ca3fa75bSLois Curfman McInnes #undef __FUNC__ 1475615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_MPIDense" 1488f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array) 149ff14e315SSatish Balay { 1503a40ed3dSBarry Smith PetscFunctionBegin; 1513a40ed3dSBarry Smith PetscFunctionReturn(0); 152ff14e315SSatish Balay } 153ff14e315SSatish Balay 1545615d1e5SSatish Balay #undef __FUNC__ 1555615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIDense" 1568f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 1578965ea79SLois Curfman McInnes { 15839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1598965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 160d36fbae8SSatish Balay int ierr,nstash,reallocs; 1618965ea79SLois Curfman McInnes InsertMode addv; 1628965ea79SLois Curfman McInnes 1633a40ed3dSBarry Smith PetscFunctionBegin; 1648965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 165ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1667056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 167a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Cannot mix adds/inserts on different procs"); 1688965ea79SLois Curfman McInnes } 169e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1708965ea79SLois Curfman McInnes 1718798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 1728798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 173d36fbae8SSatish Balay PLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n", 174d36fbae8SSatish Balay nstash,reallocs); 1753a40ed3dSBarry Smith PetscFunctionReturn(0); 1768965ea79SLois Curfman McInnes } 1778965ea79SLois Curfman McInnes 1785615d1e5SSatish Balay #undef __FUNC__ 1795615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIDense" 1808f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 1818965ea79SLois Curfman McInnes { 18239ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 1837ef1d9bdSSatish Balay int i,n,ierr,*row,*col,flg,j,rstart,ncols; 1847ef1d9bdSSatish Balay Scalar *val; 185e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 1868965ea79SLois Curfman McInnes 1873a40ed3dSBarry Smith PetscFunctionBegin; 1888965ea79SLois Curfman McInnes /* wait on receives */ 1897ef1d9bdSSatish Balay while (1) { 1908798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 1917ef1d9bdSSatish Balay if (!flg) break; 1928965ea79SLois Curfman McInnes 1937ef1d9bdSSatish Balay for ( i=0; i<n; ) { 1947ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 1957ef1d9bdSSatish Balay for ( j=i,rstart=row[j]; j<n; j++ ) { if (row[j] != rstart) break; } 1967ef1d9bdSSatish Balay if (j < n) ncols = j-i; 1977ef1d9bdSSatish Balay else ncols = n-i; 1987ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 1997ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2007ef1d9bdSSatish Balay i = j; 2018965ea79SLois Curfman McInnes } 2027ef1d9bdSSatish Balay } 2038798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2048965ea79SLois Curfman McInnes 20539ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 20639ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2078965ea79SLois Curfman McInnes 2086d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 20939ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2108965ea79SLois Curfman McInnes } 2113a40ed3dSBarry Smith PetscFunctionReturn(0); 2128965ea79SLois Curfman McInnes } 2138965ea79SLois Curfman McInnes 2145615d1e5SSatish Balay #undef __FUNC__ 2155615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIDense" 2168f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A) 2178965ea79SLois Curfman McInnes { 2183a40ed3dSBarry Smith int ierr; 21939ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 2203a40ed3dSBarry Smith 2213a40ed3dSBarry Smith PetscFunctionBegin; 2223a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2233a40ed3dSBarry Smith PetscFunctionReturn(0); 2248965ea79SLois Curfman McInnes } 2258965ea79SLois Curfman McInnes 2265615d1e5SSatish Balay #undef __FUNC__ 2275615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_MPIDense" 2288f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs) 2294e220ebcSLois Curfman McInnes { 2303a40ed3dSBarry Smith PetscFunctionBegin; 2314e220ebcSLois Curfman McInnes *bs = 1; 2323a40ed3dSBarry Smith PetscFunctionReturn(0); 2334e220ebcSLois Curfman McInnes } 2344e220ebcSLois Curfman McInnes 2358965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2368965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2378965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2383501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2398965ea79SLois Curfman McInnes routine. 2408965ea79SLois Curfman McInnes */ 2415615d1e5SSatish Balay #undef __FUNC__ 2425615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIDense" 2438f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 2448965ea79SLois Curfman McInnes { 24539ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 2468965ea79SLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 2478965ea79SLois Curfman McInnes int *procs,*nprocs,j,found,idx,nsends,*work; 2488965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2498965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2508965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2518965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2528965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2538965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2548965ea79SLois Curfman McInnes IS istmp; 2558965ea79SLois Curfman McInnes 2563a40ed3dSBarry Smith PetscFunctionBegin; 25777c4ece6SBarry Smith ierr = ISGetSize(is,&N);CHKERRQ(ierr); 2588965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 2598965ea79SLois Curfman McInnes 2608965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 2610452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) );CHKPTRQ(nprocs); 262549d3d68SSatish Balay ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 263549d3d68SSatish Balay procs = nprocs + size; 2640452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int));CHKPTRQ(owner); /* see note*/ 2658965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2668965ea79SLois Curfman McInnes idx = rows[i]; 2678965ea79SLois Curfman McInnes found = 0; 2688965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 2698965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 2708965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 2718965ea79SLois Curfman McInnes } 2728965ea79SLois Curfman McInnes } 273a8c6a408SBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range"); 2748965ea79SLois Curfman McInnes } 2758965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 2768965ea79SLois Curfman McInnes 2778965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 2780452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) );CHKPTRQ(work); 279ca161407SBarry Smith ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 2808965ea79SLois Curfman McInnes nrecvs = work[rank]; 281ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 2828965ea79SLois Curfman McInnes nmax = work[rank]; 283606d414cSSatish Balay ierr = PetscFree(work);CHKERRQ(ierr); 2848965ea79SLois Curfman McInnes 2858965ea79SLois Curfman McInnes /* post receives: */ 2863a40ed3dSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 2873a40ed3dSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 2888965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 289ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 2908965ea79SLois Curfman McInnes } 2918965ea79SLois Curfman McInnes 2928965ea79SLois Curfman McInnes /* do sends: 2938965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 2948965ea79SLois Curfman McInnes the ith processor 2958965ea79SLois Curfman McInnes */ 2960452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) );CHKPTRQ(svalues); 2977056b6fcSBarry Smith send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 2980452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) );CHKPTRQ(starts); 2998965ea79SLois Curfman McInnes starts[0] = 0; 3008965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3018965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 3028965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3038965ea79SLois Curfman McInnes } 3048965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3058965ea79SLois Curfman McInnes 3068965ea79SLois Curfman McInnes starts[0] = 0; 3078965ea79SLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3088965ea79SLois Curfman McInnes count = 0; 3098965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 3108965ea79SLois Curfman McInnes if (procs[i]) { 311ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3128965ea79SLois Curfman McInnes } 3138965ea79SLois Curfman McInnes } 314606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3158965ea79SLois Curfman McInnes 3168965ea79SLois Curfman McInnes base = owners[rank]; 3178965ea79SLois Curfman McInnes 3188965ea79SLois Curfman McInnes /* wait on receives */ 3190452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) );CHKPTRQ(lens); 3208965ea79SLois Curfman McInnes source = lens + nrecvs; 3218965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3228965ea79SLois Curfman McInnes while (count) { 323ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3248965ea79SLois Curfman McInnes /* unpack receives into our local space */ 325ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 3268965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3278965ea79SLois Curfman McInnes lens[imdex] = n; 3288965ea79SLois Curfman McInnes slen += n; 3298965ea79SLois Curfman McInnes count--; 3308965ea79SLois Curfman McInnes } 331606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3328965ea79SLois Curfman McInnes 3338965ea79SLois Curfman McInnes /* move the data into the send scatter */ 3340452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) );CHKPTRQ(lrows); 3358965ea79SLois Curfman McInnes count = 0; 3368965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 3378965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3388965ea79SLois Curfman McInnes for ( j=0; j<lens[i]; j++ ) { 3398965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3408965ea79SLois Curfman McInnes } 3418965ea79SLois Curfman McInnes } 342606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 343606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 344606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 345606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3468965ea79SLois Curfman McInnes 3478965ea79SLois Curfman McInnes /* actually zap the local rows */ 348029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 3498965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 350606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3518965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 3528965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 3538965ea79SLois Curfman McInnes 3548965ea79SLois Curfman McInnes /* wait on sends */ 3558965ea79SLois Curfman McInnes if (nsends) { 3567056b6fcSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 357ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 358606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 3598965ea79SLois Curfman McInnes } 360606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 361606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 3628965ea79SLois Curfman McInnes 3633a40ed3dSBarry Smith PetscFunctionReturn(0); 3648965ea79SLois Curfman McInnes } 3658965ea79SLois Curfman McInnes 3665615d1e5SSatish Balay #undef __FUNC__ 3675615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIDense" 3688f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3698965ea79SLois Curfman McInnes { 37039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3718965ea79SLois Curfman McInnes int ierr; 372c456f294SBarry Smith 3733a40ed3dSBarry Smith PetscFunctionBegin; 37443a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 37543a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 37644cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 3773a40ed3dSBarry Smith PetscFunctionReturn(0); 3788965ea79SLois Curfman McInnes } 3798965ea79SLois Curfman McInnes 3805615d1e5SSatish Balay #undef __FUNC__ 3815615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIDense" 3828f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3838965ea79SLois Curfman McInnes { 38439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3858965ea79SLois Curfman McInnes int ierr; 386c456f294SBarry Smith 3873a40ed3dSBarry Smith PetscFunctionBegin; 38843a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 38943a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39044cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 3913a40ed3dSBarry Smith PetscFunctionReturn(0); 3928965ea79SLois Curfman McInnes } 3938965ea79SLois Curfman McInnes 3945615d1e5SSatish Balay #undef __FUNC__ 3955615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIDense" 3968f6be9afSLois Curfman McInnes int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 397096963f5SLois Curfman McInnes { 398096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 399096963f5SLois Curfman McInnes int ierr; 4003501a2bdSLois Curfman McInnes Scalar zero = 0.0; 401096963f5SLois Curfman McInnes 4023a40ed3dSBarry Smith PetscFunctionBegin; 4033501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 40444cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 405537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 406537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4073a40ed3dSBarry Smith PetscFunctionReturn(0); 408096963f5SLois Curfman McInnes } 409096963f5SLois Curfman McInnes 4105615d1e5SSatish Balay #undef __FUNC__ 4115615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIDense" 4128f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 413096963f5SLois Curfman McInnes { 414096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 415096963f5SLois Curfman McInnes int ierr; 416096963f5SLois Curfman McInnes 4173a40ed3dSBarry Smith PetscFunctionBegin; 4183501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 41944cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 420537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 421537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4223a40ed3dSBarry Smith PetscFunctionReturn(0); 423096963f5SLois Curfman McInnes } 424096963f5SLois Curfman McInnes 4255615d1e5SSatish Balay #undef __FUNC__ 4265615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIDense" 4278f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 4288965ea79SLois Curfman McInnes { 42939ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 430096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 43144cd7ae7SLois Curfman McInnes int ierr, len, i, n, m = a->m, radd; 43244cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 433ed3cc1f0SBarry Smith 4343a40ed3dSBarry Smith PetscFunctionBegin; 43544cd7ae7SLois Curfman McInnes VecSet(&zero,v); 436096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x);CHKERRQ(ierr); 437096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 438a8c6a408SBarry Smith if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 43944cd7ae7SLois Curfman McInnes len = PetscMin(aloc->m,aloc->n); 4407ddc982cSLois Curfman McInnes radd = a->rstart*m; 44144cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 442096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 443096963f5SLois Curfman McInnes } 4449a8c540fSSatish Balay ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 4453a40ed3dSBarry Smith PetscFunctionReturn(0); 4468965ea79SLois Curfman McInnes } 4478965ea79SLois Curfman McInnes 4485615d1e5SSatish Balay #undef __FUNC__ 4495615d1e5SSatish Balay #define __FUNC__ "MatDestroy_MPIDense" 450e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 4518965ea79SLois Curfman McInnes { 4523501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4538965ea79SLois Curfman McInnes int ierr; 454ed3cc1f0SBarry Smith 4553a40ed3dSBarry Smith PetscFunctionBegin; 45694d884c6SBarry Smith 45794d884c6SBarry Smith if (mat->mapping) { 45894d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 45994d884c6SBarry Smith } 46094d884c6SBarry Smith if (mat->bmapping) { 46194d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 46294d884c6SBarry Smith } 463aa482453SBarry Smith #if defined(PETSC_USE_LOG) 464e1311b90SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N); 4658965ea79SLois Curfman McInnes #endif 4668798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 467606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 4683501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 4693501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4703501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 471622d7880SLois Curfman McInnes if (mdn->factor) { 472606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 473606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 474606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 475606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 476622d7880SLois Curfman McInnes } 477606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 47861b13de0SBarry Smith if (mat->rmap) { 47961b13de0SBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 48061b13de0SBarry Smith } 48161b13de0SBarry Smith if (mat->cmap) { 48261b13de0SBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 48390f02eecSBarry Smith } 4848965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 4850452661fSBarry Smith PetscHeaderDestroy(mat); 4863a40ed3dSBarry Smith PetscFunctionReturn(0); 4878965ea79SLois Curfman McInnes } 48839ddd567SLois Curfman McInnes 4895615d1e5SSatish Balay #undef __FUNC__ 4905615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary" 49139ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 4928965ea79SLois Curfman McInnes { 49339ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4948965ea79SLois Curfman McInnes int ierr; 4957056b6fcSBarry Smith 4963a40ed3dSBarry Smith PetscFunctionBegin; 49739ddd567SLois Curfman McInnes if (mdn->size == 1) { 49839ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 4998965ea79SLois Curfman McInnes } 500a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 5013a40ed3dSBarry Smith PetscFunctionReturn(0); 5028965ea79SLois Curfman McInnes } 5038965ea79SLois Curfman McInnes 5045615d1e5SSatish Balay #undef __FUNC__ 5055615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_ASCII" 50639ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 5078965ea79SLois Curfman McInnes { 50839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 50977ed5343SBarry Smith int ierr, format, size = mdn->size, rank = mdn->rank; 5108965ea79SLois Curfman McInnes FILE *fd; 51119bcc07fSBarry Smith ViewerType vtype; 5128965ea79SLois Curfman McInnes 5133a40ed3dSBarry Smith PetscFunctionBegin; 5143a40ed3dSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 51590ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 516888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 517639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 5184e220ebcSLois Curfman McInnes MatInfo info; 519888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 52077c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 5214e220ebcSLois Curfman McInnes fprintf(fd," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m, 5224e220ebcSLois Curfman McInnes (int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 523096963f5SLois Curfman McInnes fflush(fd); 52477c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 5253501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5263a40ed3dSBarry Smith PetscFunctionReturn(0); 52796f6c058SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO) { 5283a40ed3dSBarry Smith PetscFunctionReturn(0); 5298965ea79SLois Curfman McInnes } 53077ed5343SBarry Smith 5318965ea79SLois Curfman McInnes if (size == 1) { 53239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5333a40ed3dSBarry Smith } else { 5348965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5358965ea79SLois Curfman McInnes Mat A; 53639ddd567SLois Curfman McInnes int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 53739ddd567SLois Curfman McInnes Scalar *vals; 53839ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 5398965ea79SLois Curfman McInnes 5408965ea79SLois Curfman McInnes if (!rank) { 5410513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 5423a40ed3dSBarry Smith } else { 5430513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,0,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 5448965ea79SLois Curfman McInnes } 5458965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 5468965ea79SLois Curfman McInnes 54739ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 54839ddd567SLois Curfman McInnes but it's quick for now */ 54939ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 5508965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 55139ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 55239ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 55339ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 55439ddd567SLois Curfman McInnes row++; 5558965ea79SLois Curfman McInnes } 5568965ea79SLois Curfman McInnes 5576d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5586d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5598965ea79SLois Curfman McInnes if (!rank) { 56039ddd567SLois Curfman McInnes ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer);CHKERRQ(ierr); 5618965ea79SLois Curfman McInnes } 5628965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 5638965ea79SLois Curfman McInnes } 5643a40ed3dSBarry Smith PetscFunctionReturn(0); 5658965ea79SLois Curfman McInnes } 5668965ea79SLois Curfman McInnes 5675615d1e5SSatish Balay #undef __FUNC__ 5685615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense" 569e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer) 5708965ea79SLois Curfman McInnes { 57139ddd567SLois Curfman McInnes int ierr; 572bcd2baecSBarry Smith ViewerType vtype; 5738965ea79SLois Curfman McInnes 574bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 5753f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 57639ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer);CHKERRQ(ierr); 5773f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 5783a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 5795cd90555SBarry Smith } else { 5805cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 5818965ea79SLois Curfman McInnes } 5823a40ed3dSBarry Smith PetscFunctionReturn(0); 5838965ea79SLois Curfman McInnes } 5848965ea79SLois Curfman McInnes 5855615d1e5SSatish Balay #undef __FUNC__ 5865615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense" 5878f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 5888965ea79SLois Curfman McInnes { 5893501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 5903501a2bdSLois Curfman McInnes Mat mdn = mat->A; 5914e220ebcSLois Curfman McInnes int ierr; 5924e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 5938965ea79SLois Curfman McInnes 5943a40ed3dSBarry Smith PetscFunctionBegin; 5954e220ebcSLois Curfman McInnes info->rows_global = (double)mat->M; 5964e220ebcSLois Curfman McInnes info->columns_global = (double)mat->N; 5974e220ebcSLois Curfman McInnes info->rows_local = (double)mat->m; 5984e220ebcSLois Curfman McInnes info->columns_local = (double)mat->N; 5994e220ebcSLois Curfman McInnes info->block_size = 1.0; 6004e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6014e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6024e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6038965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6044e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6054e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6064e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6074e220ebcSLois Curfman McInnes info->memory = isend[3]; 6084e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6098965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 610f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 6114e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6124e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6134e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6144e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6154e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6168965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 617f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 6184e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6194e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6204e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6214e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6224e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6238965ea79SLois Curfman McInnes } 6244e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6254e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6264e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6273a40ed3dSBarry Smith PetscFunctionReturn(0); 6288965ea79SLois Curfman McInnes } 6298965ea79SLois Curfman McInnes 6308c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*); 6318aaee692SLois Curfman McInnes extern int MatLUFactorNumeric_MPIDense(Mat,Mat*); 6328aaee692SLois Curfman McInnes extern int MatLUFactor_MPIDense(Mat,IS,IS,double); 6338aaee692SLois Curfman McInnes extern int MatSolve_MPIDense(Mat,Vec,Vec); 6348c469469SLois Curfman McInnes extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec); 6358aaee692SLois Curfman McInnes extern int MatSolveTrans_MPIDense(Mat,Vec,Vec); 6368aaee692SLois Curfman McInnes extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */ 6378aaee692SLois Curfman McInnes 6385615d1e5SSatish Balay #undef __FUNC__ 6395615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense" 6408f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 6418965ea79SLois Curfman McInnes { 64239ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 6438965ea79SLois Curfman McInnes 6443a40ed3dSBarry Smith PetscFunctionBegin; 6456d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 6466d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 6474787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 6484787f768SSatish Balay op == MAT_NEW_NONZERO_ALLOCATION_ERR || 649219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 650219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED) { 651b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 652b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 653aeafbbfcSLois Curfman McInnes a->roworiented = 1; 6548965ea79SLois Curfman McInnes MatSetOption(a->A,op); 655b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 656219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 6576d4a8577SBarry Smith op == MAT_SYMMETRIC || 6586d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 659b51ba29fSSatish Balay op == MAT_YES_NEW_DIAGONALS || 660b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 661981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 6623a40ed3dSBarry Smith } else if (op == MAT_COLUMN_ORIENTED) { 6633a40ed3dSBarry Smith a->roworiented = 0; MatSetOption(a->A,op); 6643782ba37SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 6653782ba37SSatish Balay a->donotstash = 1; 6663a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS) { 6673a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 6683a40ed3dSBarry Smith } else { 6693a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 6703a40ed3dSBarry Smith } 6713a40ed3dSBarry Smith PetscFunctionReturn(0); 6728965ea79SLois Curfman McInnes } 6738965ea79SLois Curfman McInnes 6745615d1e5SSatish Balay #undef __FUNC__ 6755615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense" 6768f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n) 6778965ea79SLois Curfman McInnes { 6783501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6793a40ed3dSBarry Smith 6803a40ed3dSBarry Smith PetscFunctionBegin; 6818965ea79SLois Curfman McInnes *m = mat->M; *n = mat->N; 6823a40ed3dSBarry Smith PetscFunctionReturn(0); 6838965ea79SLois Curfman McInnes } 6848965ea79SLois Curfman McInnes 6855615d1e5SSatish Balay #undef __FUNC__ 6865615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense" 6878f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 6888965ea79SLois Curfman McInnes { 6893501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6903a40ed3dSBarry Smith 6913a40ed3dSBarry Smith PetscFunctionBegin; 6928965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 6933a40ed3dSBarry Smith PetscFunctionReturn(0); 6948965ea79SLois Curfman McInnes } 6958965ea79SLois Curfman McInnes 6965615d1e5SSatish Balay #undef __FUNC__ 6975615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense" 6988f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 6998965ea79SLois Curfman McInnes { 7003501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 7013a40ed3dSBarry Smith 7023a40ed3dSBarry Smith PetscFunctionBegin; 7038965ea79SLois Curfman McInnes *m = mat->rstart; *n = mat->rend; 7043a40ed3dSBarry Smith PetscFunctionReturn(0); 7058965ea79SLois Curfman McInnes } 7068965ea79SLois Curfman McInnes 7075615d1e5SSatish Balay #undef __FUNC__ 7085615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense" 7098f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 7108965ea79SLois Curfman McInnes { 7113501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 7123a40ed3dSBarry Smith int lrow, rstart = mat->rstart, rend = mat->rend,ierr; 7138965ea79SLois Curfman McInnes 7143a40ed3dSBarry Smith PetscFunctionBegin; 715a8c6a408SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows") 7168965ea79SLois Curfman McInnes lrow = row - rstart; 7173a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 7183a40ed3dSBarry Smith PetscFunctionReturn(0); 7198965ea79SLois Curfman McInnes } 7208965ea79SLois Curfman McInnes 7215615d1e5SSatish Balay #undef __FUNC__ 7225615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense" 7238f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 7248965ea79SLois Curfman McInnes { 725606d414cSSatish Balay int ierr; 726606d414cSSatish Balay 7273a40ed3dSBarry Smith PetscFunctionBegin; 728606d414cSSatish Balay if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 729606d414cSSatish Balay if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 7303a40ed3dSBarry Smith PetscFunctionReturn(0); 7318965ea79SLois Curfman McInnes } 7328965ea79SLois Curfman McInnes 7335615d1e5SSatish Balay #undef __FUNC__ 7345615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense" 7358f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm) 736096963f5SLois Curfman McInnes { 7373501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) A->data; 7383501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data; 7393501a2bdSLois Curfman McInnes int ierr, i, j; 7403501a2bdSLois Curfman McInnes double sum = 0.0; 7413501a2bdSLois Curfman McInnes Scalar *v = mat->v; 7423501a2bdSLois Curfman McInnes 7433a40ed3dSBarry Smith PetscFunctionBegin; 7443501a2bdSLois Curfman McInnes if (mdn->size == 1) { 7453501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm);CHKERRQ(ierr); 7463501a2bdSLois Curfman McInnes } else { 7473501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 7483501a2bdSLois Curfman McInnes for (i=0; i<mat->n*mat->m; i++ ) { 749aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 750e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 7513501a2bdSLois Curfman McInnes #else 7523501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 7533501a2bdSLois Curfman McInnes #endif 7543501a2bdSLois Curfman McInnes } 755ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 7563501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 7573501a2bdSLois Curfman McInnes PLogFlops(2*mat->n*mat->m); 7583a40ed3dSBarry Smith } else if (type == NORM_1) { 7593501a2bdSLois Curfman McInnes double *tmp, *tmp2; 7600452661fSBarry Smith tmp = (double *) PetscMalloc( 2*mdn->N*sizeof(double) );CHKPTRQ(tmp); 7613501a2bdSLois Curfman McInnes tmp2 = tmp + mdn->N; 762549d3d68SSatish Balay ierr = PetscMemzero(tmp,2*mdn->N*sizeof(double));CHKERRQ(ierr); 763096963f5SLois Curfman McInnes *norm = 0.0; 7643501a2bdSLois Curfman McInnes v = mat->v; 7653501a2bdSLois Curfman McInnes for ( j=0; j<mat->n; j++ ) { 7663501a2bdSLois Curfman McInnes for ( i=0; i<mat->m; i++ ) { 76767e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 7683501a2bdSLois Curfman McInnes } 7693501a2bdSLois Curfman McInnes } 770ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 7713501a2bdSLois Curfman McInnes for ( j=0; j<mdn->N; j++ ) { 7723501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 7733501a2bdSLois Curfman McInnes } 774606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 7753501a2bdSLois Curfman McInnes PLogFlops(mat->n*mat->m); 7763a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 7773501a2bdSLois Curfman McInnes double ntemp; 7783501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 779ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 7803a40ed3dSBarry Smith } else { 781a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 7823501a2bdSLois Curfman McInnes } 7833501a2bdSLois Curfman McInnes } 7843a40ed3dSBarry Smith PetscFunctionReturn(0); 7853501a2bdSLois Curfman McInnes } 7863501a2bdSLois Curfman McInnes 7875615d1e5SSatish Balay #undef __FUNC__ 7885615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense" 7898f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 7903501a2bdSLois Curfman McInnes { 7913501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 7923501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data; 7933501a2bdSLois Curfman McInnes Mat B; 7943501a2bdSLois Curfman McInnes int M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart; 7953501a2bdSLois Curfman McInnes int j, i, ierr; 7963501a2bdSLois Curfman McInnes Scalar *v; 7973501a2bdSLois Curfman McInnes 7983a40ed3dSBarry Smith PetscFunctionBegin; 7997056b6fcSBarry Smith if (matout == PETSC_NULL && M != N) { 800a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place"); 8017056b6fcSBarry Smith } 8027056b6fcSBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 8033501a2bdSLois Curfman McInnes 8043501a2bdSLois Curfman McInnes m = Aloc->m; n = Aloc->n; v = Aloc->v; 8050452661fSBarry Smith rwork = (int *) PetscMalloc(n*sizeof(int));CHKPTRQ(rwork); 8063501a2bdSLois Curfman McInnes for ( j=0; j<n; j++ ) { 8073501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8083501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8093501a2bdSLois Curfman McInnes v += m; 8103501a2bdSLois Curfman McInnes } 811606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 8126d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8136d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8143638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 8153501a2bdSLois Curfman McInnes *matout = B; 8163501a2bdSLois Curfman McInnes } else { 817f830108cSBarry Smith PetscOps *Abops; 81809dc0095SBarry Smith MatOps Aops; 819f830108cSBarry Smith 8203501a2bdSLois Curfman McInnes /* This isn't really an in-place transpose, but free data struct from a */ 821606d414cSSatish Balay ierr = PetscFree(a->rowners);CHKERRQ(ierr); 8223501a2bdSLois Curfman McInnes ierr = MatDestroy(a->A);CHKERRQ(ierr); 8233501a2bdSLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 8243501a2bdSLois Curfman McInnes if (a->Mvctx) VecScatterDestroy(a->Mvctx); 825606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 826f830108cSBarry Smith 827f830108cSBarry Smith /* 828f830108cSBarry Smith This is horrible, horrible code. We need to keep the 829f830108cSBarry Smith A pointers for the bops and ops but copy everything 830f830108cSBarry Smith else from C. 831f830108cSBarry Smith */ 832f830108cSBarry Smith Abops = A->bops; 833f830108cSBarry Smith Aops = A->ops; 834549d3d68SSatish Balay ierr = PetscMemcpy(A,B,sizeof(struct _p_Mat));CHKERRQ(ierr); 835f830108cSBarry Smith A->bops = Abops; 836f830108cSBarry Smith A->ops = Aops; 837f830108cSBarry Smith 8380452661fSBarry Smith PetscHeaderDestroy(B); 8393501a2bdSLois Curfman McInnes } 8403a40ed3dSBarry Smith PetscFunctionReturn(0); 841096963f5SLois Curfman McInnes } 842096963f5SLois Curfman McInnes 843eadb2fb4SBarry Smith #include "pinclude/blaslapack.h" 8445615d1e5SSatish Balay #undef __FUNC__ 8455615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense" 8468f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA) 84744cd7ae7SLois Curfman McInnes { 84844cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense *) inA->data; 84944cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->A->data; 85044cd7ae7SLois Curfman McInnes int one = 1, nz; 85144cd7ae7SLois Curfman McInnes 8523a40ed3dSBarry Smith PetscFunctionBegin; 85344cd7ae7SLois Curfman McInnes nz = a->m*a->n; 85444cd7ae7SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 85544cd7ae7SLois Curfman McInnes PLogFlops(nz); 8563a40ed3dSBarry Smith PetscFunctionReturn(0); 85744cd7ae7SLois Curfman McInnes } 85844cd7ae7SLois Curfman McInnes 8595609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 8607b2a1423SBarry Smith extern int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **); 8618965ea79SLois Curfman McInnes 8628965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 86309dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 86409dc0095SBarry Smith MatGetRow_MPIDense, 86509dc0095SBarry Smith MatRestoreRow_MPIDense, 86609dc0095SBarry Smith MatMult_MPIDense, 86709dc0095SBarry Smith MatMultAdd_MPIDense, 86809dc0095SBarry Smith MatMultTrans_MPIDense, 86909dc0095SBarry Smith MatMultTransAdd_MPIDense, 8708965ea79SLois Curfman McInnes 0, 87109dc0095SBarry Smith 0, 87209dc0095SBarry Smith 0, 87309dc0095SBarry Smith 0, 87409dc0095SBarry Smith 0, 87509dc0095SBarry Smith 0, 87609dc0095SBarry Smith 0, 87709dc0095SBarry Smith MatTranspose_MPIDense, 87809dc0095SBarry Smith MatGetInfo_MPIDense,0, 87909dc0095SBarry Smith MatGetDiagonal_MPIDense, 88009dc0095SBarry Smith 0, 88109dc0095SBarry Smith MatNorm_MPIDense, 88209dc0095SBarry Smith MatAssemblyBegin_MPIDense, 88309dc0095SBarry Smith MatAssemblyEnd_MPIDense, 88409dc0095SBarry Smith 0, 88509dc0095SBarry Smith MatSetOption_MPIDense, 88609dc0095SBarry Smith MatZeroEntries_MPIDense, 88709dc0095SBarry Smith MatZeroRows_MPIDense, 88809dc0095SBarry Smith 0, 88909dc0095SBarry Smith 0, 89009dc0095SBarry Smith 0, 89109dc0095SBarry Smith 0, 89209dc0095SBarry Smith MatGetSize_MPIDense, 89309dc0095SBarry Smith MatGetLocalSize_MPIDense, 89439ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 89509dc0095SBarry Smith 0, 89609dc0095SBarry Smith 0, 89709dc0095SBarry Smith MatGetArray_MPIDense, 89809dc0095SBarry Smith MatRestoreArray_MPIDense, 8995609ef8eSBarry Smith MatDuplicate_MPIDense, 90009dc0095SBarry Smith 0, 90109dc0095SBarry Smith 0, 90209dc0095SBarry Smith 0, 90309dc0095SBarry Smith 0, 90409dc0095SBarry Smith 0, 9052ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 90609dc0095SBarry Smith 0, 90709dc0095SBarry Smith MatGetValues_MPIDense, 90809dc0095SBarry Smith 0, 90909dc0095SBarry Smith 0, 91009dc0095SBarry Smith MatScale_MPIDense, 91109dc0095SBarry Smith 0, 91209dc0095SBarry Smith 0, 91309dc0095SBarry Smith 0, 91409dc0095SBarry Smith MatGetBlockSize_MPIDense, 91509dc0095SBarry Smith 0, 91609dc0095SBarry Smith 0, 91709dc0095SBarry Smith 0, 91809dc0095SBarry Smith 0, 91909dc0095SBarry Smith 0, 92009dc0095SBarry Smith 0, 92109dc0095SBarry Smith 0, 92209dc0095SBarry Smith 0, 92309dc0095SBarry Smith 0, 924ca3fa75bSLois Curfman McInnes MatGetSubMatrix_MPIDense, 92509dc0095SBarry Smith 0, 92609dc0095SBarry Smith 0, 92709dc0095SBarry Smith MatGetMaps_Petsc}; 9288965ea79SLois Curfman McInnes 9295615d1e5SSatish Balay #undef __FUNC__ 9305615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense" 9318965ea79SLois Curfman McInnes /*@C 93239ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 9338965ea79SLois Curfman McInnes 934db81eaa0SLois Curfman McInnes Collective on MPI_Comm 935db81eaa0SLois Curfman McInnes 9368965ea79SLois Curfman McInnes Input Parameters: 937db81eaa0SLois Curfman McInnes + comm - MPI communicator 9388965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 939db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 9408965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 941db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 942db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 943dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 9448965ea79SLois Curfman McInnes 9458965ea79SLois Curfman McInnes Output Parameter: 946477f1c0bSLois Curfman McInnes . A - the matrix 9478965ea79SLois Curfman McInnes 948b259b22eSLois Curfman McInnes Notes: 94939ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 95039ddd567SLois Curfman McInnes storage by columns. 9518965ea79SLois Curfman McInnes 95218f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 95318f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 954b4fd4287SBarry Smith set data=PETSC_NULL. 95518f449edSLois Curfman McInnes 9568965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 9578965ea79SLois Curfman McInnes (possibly both). 9588965ea79SLois Curfman McInnes 9593501a2bdSLois Curfman McInnes Currently, the only parallel dense matrix decomposition is by rows, 9603501a2bdSLois Curfman McInnes so that n=N and each submatrix owns all of the global columns. 9613501a2bdSLois Curfman McInnes 962027ccd11SLois Curfman McInnes Level: intermediate 963027ccd11SLois Curfman McInnes 96439ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel 9658965ea79SLois Curfman McInnes 96639ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 9678965ea79SLois Curfman McInnes @*/ 968477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A) 9698965ea79SLois Curfman McInnes { 9708965ea79SLois Curfman McInnes Mat mat; 97139ddd567SLois Curfman McInnes Mat_MPIDense *a; 97225cdf11fSBarry Smith int ierr, i,flg; 9738965ea79SLois Curfman McInnes 9743a40ed3dSBarry Smith PetscFunctionBegin; 975ed2daf61SLois Curfman McInnes /* Note: For now, when data is specified above, this assumes the user correctly 976ed2daf61SLois Curfman McInnes allocates the local dense storage space. We should add error checking. */ 97718f449edSLois Curfman McInnes 978477f1c0bSLois Curfman McInnes *A = 0; 9793f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",comm,MatDestroy,MatView); 9808965ea79SLois Curfman McInnes PLogObjectCreate(mat); 9810452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 982549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 983e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 984e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 9858965ea79SLois Curfman McInnes mat->factor = 0; 98690f02eecSBarry Smith mat->mapping = 0; 9878965ea79SLois Curfman McInnes 988622d7880SLois Curfman McInnes a->factor = 0; 989e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 990d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&a->rank);CHKERRQ(ierr); 991d132466eSBarry Smith ierr = MPI_Comm_size(comm,&a->size);CHKERRQ(ierr); 9928965ea79SLois Curfman McInnes 99396f6c058SBarry Smith ierr = PetscSplitOwnership(comm,&m,&M);CHKERRQ(ierr); 99439ddd567SLois Curfman McInnes 995c7fcc2eaSBarry Smith /* 996c7fcc2eaSBarry Smith The computation of n is wrong below, n should represent the number of local 997c7fcc2eaSBarry Smith rows in the right (column vector) 998c7fcc2eaSBarry Smith */ 999c7fcc2eaSBarry Smith 100039ddd567SLois Curfman McInnes /* each row stores all columns */ 100139ddd567SLois Curfman McInnes if (N == PETSC_DECIDE) N = n; 1002d7e8b826SBarry Smith if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);} 1003a8c6a408SBarry Smith /* if (n != N) SETERRQ(PETSC_ERR_SUP,0,"For now, only n=N is supported"); */ 1004aca0ad90SLois Curfman McInnes a->N = mat->N = N; 1005aca0ad90SLois Curfman McInnes a->M = mat->M = M; 1006aca0ad90SLois Curfman McInnes a->m = mat->m = m; 1007aca0ad90SLois Curfman McInnes a->n = mat->n = n; 10088965ea79SLois Curfman McInnes 1009c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 1010c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 1011488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,M,&mat->rmap);CHKERRQ(ierr); 101296f6c058SBarry Smith ierr = MapCreateMPI(comm,PETSC_DECIDE,N,&mat->cmap);CHKERRQ(ierr); 1013c7fcc2eaSBarry Smith 10148965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 1015d7e8b826SBarry Smith a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners); 1016d7e8b826SBarry Smith a->cowners = a->rowners + a->size + 1; 1017f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1018ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 10198965ea79SLois Curfman McInnes a->rowners[0] = 0; 10208965ea79SLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 10218965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 10228965ea79SLois Curfman McInnes } 10238965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 10248965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 1025ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1026d7e8b826SBarry Smith a->cowners[0] = 0; 1027d7e8b826SBarry Smith for ( i=2; i<=a->size; i++ ) { 1028d7e8b826SBarry Smith a->cowners[i] += a->cowners[i-1]; 1029d7e8b826SBarry Smith } 10308965ea79SLois Curfman McInnes 1031029af93fSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A);CHKERRQ(ierr); 10328965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 10338965ea79SLois Curfman McInnes 10348965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 10353782ba37SSatish Balay a->donotstash = 0; 10368798bf22SSatish Balay ierr = MatStashCreate_Private(comm,1,&mat->stash);CHKERRQ(ierr); 10378965ea79SLois Curfman McInnes 10388965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 10398965ea79SLois Curfman McInnes a->lvec = 0; 10408965ea79SLois Curfman McInnes a->Mvctx = 0; 104139b7565bSBarry Smith a->roworiented = 1; 10428965ea79SLois Curfman McInnes 1043477f1c0bSLois Curfman McInnes *A = mat; 104425cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 104525cdf11fSBarry Smith if (flg) { 10468c469469SLois Curfman McInnes ierr = MatPrintHelp(mat);CHKERRQ(ierr); 10478c469469SLois Curfman McInnes } 10483a40ed3dSBarry Smith PetscFunctionReturn(0); 10498965ea79SLois Curfman McInnes } 10508965ea79SLois Curfman McInnes 10515615d1e5SSatish Balay #undef __FUNC__ 10525609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_MPIDense" 10535609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 10548965ea79SLois Curfman McInnes { 10558965ea79SLois Curfman McInnes Mat mat; 10563501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data; 105739ddd567SLois Curfman McInnes int ierr; 10582ba99913SLois Curfman McInnes FactorCtx *factor; 10598965ea79SLois Curfman McInnes 10603a40ed3dSBarry Smith PetscFunctionBegin; 10618965ea79SLois Curfman McInnes *newmat = 0; 10623f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",A->comm,MatDestroy,MatView); 10638965ea79SLois Curfman McInnes PLogObjectCreate(mat); 10640452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 1065549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1066e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 1067e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 10683501a2bdSLois Curfman McInnes mat->factor = A->factor; 1069c456f294SBarry Smith mat->assembled = PETSC_TRUE; 10708965ea79SLois Curfman McInnes 107144cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 107244cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 107344cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 107444cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 10752ba99913SLois Curfman McInnes if (oldmat->factor) { 10762ba99913SLois Curfman McInnes a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx));CHKPTRQ(factor); 10772ba99913SLois Curfman McInnes /* copy factor contents ... add this code! */ 10782ba99913SLois Curfman McInnes } else a->factor = 0; 10798965ea79SLois Curfman McInnes 10808965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 10818965ea79SLois Curfman McInnes a->rend = oldmat->rend; 10828965ea79SLois Curfman McInnes a->size = oldmat->size; 10838965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1084e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 10853782ba37SSatish Balay a->donotstash = oldmat->donotstash; 10860452661fSBarry Smith a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners); 1087f09e8eb9SSatish Balay PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1088549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 10898798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 10908965ea79SLois Curfman McInnes 10918965ea79SLois Curfman McInnes ierr = VecDuplicate(oldmat->lvec,&a->lvec);CHKERRQ(ierr); 10928965ea79SLois Curfman McInnes PLogObjectParent(mat,a->lvec); 109355659b69SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx);CHKERRQ(ierr); 10948965ea79SLois Curfman McInnes PLogObjectParent(mat,a->Mvctx); 10955609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 10968965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 10978965ea79SLois Curfman McInnes *newmat = mat; 10983a40ed3dSBarry Smith PetscFunctionReturn(0); 10998965ea79SLois Curfman McInnes } 11008965ea79SLois Curfman McInnes 110177c4ece6SBarry Smith #include "sys.h" 11028965ea79SLois Curfman McInnes 11035615d1e5SSatish Balay #undef __FUNC__ 11045615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile" 110590ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat) 110690ace30eSBarry Smith { 110740011551SBarry Smith int *rowners, i,size,rank,m,ierr,nz,j; 110890ace30eSBarry Smith Scalar *array,*vals,*vals_ptr; 110990ace30eSBarry Smith MPI_Status status; 111090ace30eSBarry Smith 11113a40ed3dSBarry Smith PetscFunctionBegin; 1112d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1113d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 111490ace30eSBarry Smith 111590ace30eSBarry Smith /* determine ownership of all rows */ 111690ace30eSBarry Smith m = M/size + ((M % size) > rank); 111790ace30eSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1118ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 111990ace30eSBarry Smith rowners[0] = 0; 112090ace30eSBarry Smith for ( i=2; i<=size; i++ ) { 112190ace30eSBarry Smith rowners[i] += rowners[i-1]; 112290ace30eSBarry Smith } 112390ace30eSBarry Smith 112490ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 112590ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 112690ace30eSBarry Smith 112790ace30eSBarry Smith if (!rank) { 112890ace30eSBarry Smith vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals); 112990ace30eSBarry Smith 113090ace30eSBarry Smith /* read in my part of the matrix numerical values */ 11310752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 113290ace30eSBarry Smith 113390ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 113490ace30eSBarry Smith vals_ptr = vals; 113590ace30eSBarry Smith for ( i=0; i<m; i++ ) { 113690ace30eSBarry Smith for ( j=0; j<N; j++ ) { 113790ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 113890ace30eSBarry Smith } 113990ace30eSBarry Smith } 114090ace30eSBarry Smith 114190ace30eSBarry Smith /* read in other processors and ship out */ 114290ace30eSBarry Smith for ( i=1; i<size; i++ ) { 114390ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 11440752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1145ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 114690ace30eSBarry Smith } 11473a40ed3dSBarry Smith } else { 114890ace30eSBarry Smith /* receive numeric values */ 114990ace30eSBarry Smith vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals); 115090ace30eSBarry Smith 115190ace30eSBarry Smith /* receive message of values*/ 1152ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 115390ace30eSBarry Smith 115490ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 115590ace30eSBarry Smith vals_ptr = vals; 115690ace30eSBarry Smith for ( i=0; i<m; i++ ) { 115790ace30eSBarry Smith for ( j=0; j<N; j++ ) { 115890ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 115990ace30eSBarry Smith } 116090ace30eSBarry Smith } 116190ace30eSBarry Smith } 1162606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1163606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 11646d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11656d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11663a40ed3dSBarry Smith PetscFunctionReturn(0); 116790ace30eSBarry Smith } 116890ace30eSBarry Smith 116990ace30eSBarry Smith 11705615d1e5SSatish Balay #undef __FUNC__ 11715615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense" 117219bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat) 11738965ea79SLois Curfman McInnes { 11748965ea79SLois Curfman McInnes Mat A; 11758965ea79SLois Curfman McInnes Scalar *vals,*svals; 117619bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 11778965ea79SLois Curfman McInnes MPI_Status status; 11788965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 11798965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 118019bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 11813a40ed3dSBarry Smith int i, nz, ierr, j,rstart, rend, fd; 11828965ea79SLois Curfman McInnes 11833a40ed3dSBarry Smith PetscFunctionBegin; 1184d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1185d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 11868965ea79SLois Curfman McInnes if (!rank) { 118790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 11880752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1189a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 11908965ea79SLois Curfman McInnes } 11918965ea79SLois Curfman McInnes 1192ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 119390ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 119490ace30eSBarry Smith 119590ace30eSBarry Smith /* 119690ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 119790ace30eSBarry Smith */ 119890ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 11993a40ed3dSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 12003a40ed3dSBarry Smith PetscFunctionReturn(0); 120190ace30eSBarry Smith } 120290ace30eSBarry Smith 12038965ea79SLois Curfman McInnes /* determine ownership of all rows */ 12048965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 12050452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1206ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 12078965ea79SLois Curfman McInnes rowners[0] = 0; 12088965ea79SLois Curfman McInnes for ( i=2; i<=size; i++ ) { 12098965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 12108965ea79SLois Curfman McInnes } 12118965ea79SLois Curfman McInnes rstart = rowners[rank]; 12128965ea79SLois Curfman McInnes rend = rowners[rank+1]; 12138965ea79SLois Curfman McInnes 12148965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 12150452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) );CHKPTRQ(ourlens); 12168965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 12178965ea79SLois Curfman McInnes if (!rank) { 12180452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) );CHKPTRQ(rowlengths); 12190752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 12200452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(sndcounts); 12218965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1222ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1223606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1224ca161407SBarry Smith } else { 1225ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 12268965ea79SLois Curfman McInnes } 12278965ea79SLois Curfman McInnes 12288965ea79SLois Curfman McInnes if (!rank) { 12298965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 12300452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(procsnz); 1231549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 12328965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 12338965ea79SLois Curfman McInnes for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 12348965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 12358965ea79SLois Curfman McInnes } 12368965ea79SLois Curfman McInnes } 1237606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 12388965ea79SLois Curfman McInnes 12398965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 12408965ea79SLois Curfman McInnes maxnz = 0; 12418965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 12420452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 12438965ea79SLois Curfman McInnes } 12440452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) );CHKPTRQ(cols); 12458965ea79SLois Curfman McInnes 12468965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 12478965ea79SLois Curfman McInnes nz = procsnz[0]; 12480452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols); 12490752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 12508965ea79SLois Curfman McInnes 12518965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 12528965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 12538965ea79SLois Curfman McInnes nz = procsnz[i]; 12540752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1255ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 12568965ea79SLois Curfman McInnes } 1257606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 12583a40ed3dSBarry Smith } else { 12598965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 12608965ea79SLois Curfman McInnes nz = 0; 12618965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12628965ea79SLois Curfman McInnes nz += ourlens[i]; 12638965ea79SLois Curfman McInnes } 12640452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols); 12658965ea79SLois Curfman McInnes 12668965ea79SLois Curfman McInnes /* receive message of column indices*/ 1267ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1268ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1269a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12708965ea79SLois Curfman McInnes } 12718965ea79SLois Curfman McInnes 12728965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1273549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 12748965ea79SLois Curfman McInnes jj = 0; 12758965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12768965ea79SLois Curfman McInnes for ( j=0; j<ourlens[i]; j++ ) { 12778965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 12788965ea79SLois Curfman McInnes jj++; 12798965ea79SLois Curfman McInnes } 12808965ea79SLois Curfman McInnes } 12818965ea79SLois Curfman McInnes 12828965ea79SLois Curfman McInnes /* create our matrix */ 12838965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12848965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 12858965ea79SLois Curfman McInnes } 1286b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 12878965ea79SLois Curfman McInnes A = *newmat; 12888965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12898965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 12908965ea79SLois Curfman McInnes } 12918965ea79SLois Curfman McInnes 12928965ea79SLois Curfman McInnes if (!rank) { 12930452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) );CHKPTRQ(vals); 12948965ea79SLois Curfman McInnes 12958965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 12968965ea79SLois Curfman McInnes nz = procsnz[0]; 12970752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 12988965ea79SLois Curfman McInnes 12998965ea79SLois Curfman McInnes /* insert into matrix */ 13008965ea79SLois Curfman McInnes jj = rstart; 13018965ea79SLois Curfman McInnes smycols = mycols; 13028965ea79SLois Curfman McInnes svals = vals; 13038965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 13048965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 13058965ea79SLois Curfman McInnes smycols += ourlens[i]; 13068965ea79SLois Curfman McInnes svals += ourlens[i]; 13078965ea79SLois Curfman McInnes jj++; 13088965ea79SLois Curfman McInnes } 13098965ea79SLois Curfman McInnes 13108965ea79SLois Curfman McInnes /* read in other processors and ship out */ 13118965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 13128965ea79SLois Curfman McInnes nz = procsnz[i]; 13130752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1314ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 13158965ea79SLois Curfman McInnes } 1316606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 13173a40ed3dSBarry Smith } else { 13188965ea79SLois Curfman McInnes /* receive numeric values */ 13190452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) );CHKPTRQ(vals); 13208965ea79SLois Curfman McInnes 13218965ea79SLois Curfman McInnes /* receive message of values*/ 1322ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1323ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1324a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 13258965ea79SLois Curfman McInnes 13268965ea79SLois Curfman McInnes /* insert into matrix */ 13278965ea79SLois Curfman McInnes jj = rstart; 13288965ea79SLois Curfman McInnes smycols = mycols; 13298965ea79SLois Curfman McInnes svals = vals; 13308965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 13318965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 13328965ea79SLois Curfman McInnes smycols += ourlens[i]; 13338965ea79SLois Curfman McInnes svals += ourlens[i]; 13348965ea79SLois Curfman McInnes jj++; 13358965ea79SLois Curfman McInnes } 13368965ea79SLois Curfman McInnes } 1337606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1338606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1339606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1340606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 13418965ea79SLois Curfman McInnes 13426d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13436d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13443a40ed3dSBarry Smith PetscFunctionReturn(0); 13458965ea79SLois Curfman McInnes } 134690ace30eSBarry Smith 134790ace30eSBarry Smith 134890ace30eSBarry Smith 134990ace30eSBarry Smith 135090ace30eSBarry Smith 1351