1*87828ca2SBarry Smith /*$Id: mpidense.c,v 1.156 2001/07/20 21:19:32 bsmith 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 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "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; 16cfce73b9SSatish Balay int m = A->m,rstart = mdn->rstart,ierr; 17*87828ca2SBarry Smith PetscScalar *array; 180de54da6SSatish Balay MPI_Comm comm; 190de54da6SSatish Balay 200de54da6SSatish Balay PetscFunctionBegin; 21273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"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 = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 270de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 280de54da6SSatish Balay ierr = MatCreateSeqDense(comm,m,m,array+m*rstart,B);CHKERRQ(ierr); 290de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 300de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 310de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 320de54da6SSatish Balay 330de54da6SSatish Balay *iscopy = PETSC_TRUE; 340de54da6SSatish Balay PetscFunctionReturn(0); 350de54da6SSatish Balay } 360de54da6SSatish Balay EXTERN_C_END 370de54da6SSatish Balay 38ca44d042SBarry Smith EXTERN int MatSetUpMultiply_MPIDense(Mat); 397ef1d9bdSSatish Balay 404a2ae208SSatish Balay #undef __FUNCT__ 414a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 42*87828ca2SBarry Smith int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,PetscScalar *v,InsertMode addv) 438965ea79SLois Curfman McInnes { 4439b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 4539b7565bSBarry Smith int ierr,i,j,rstart = A->rstart,rend = A->rend,row; 46273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 478965ea79SLois Curfman McInnes 483a40ed3dSBarry Smith PetscFunctionBegin; 498965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 505ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 51273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 528965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 538965ea79SLois Curfman McInnes row = idxm[i] - rstart; 5439b7565bSBarry Smith if (roworiented) { 5539b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 563a40ed3dSBarry Smith } else { 578965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 585ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 59273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 6039b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 6139b7565bSBarry Smith } 628965ea79SLois Curfman McInnes } 633a40ed3dSBarry Smith } else { 643782ba37SSatish Balay if (!A->donotstash) { 6539b7565bSBarry Smith if (roworiented) { 668798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 67d36fbae8SSatish Balay } else { 688798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 6939b7565bSBarry Smith } 70b49de8d1SLois Curfman McInnes } 71b49de8d1SLois Curfman McInnes } 723782ba37SSatish Balay } 733a40ed3dSBarry Smith PetscFunctionReturn(0); 74b49de8d1SLois Curfman McInnes } 75b49de8d1SLois Curfman McInnes 764a2ae208SSatish Balay #undef __FUNCT__ 774a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 78*87828ca2SBarry Smith int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,PetscScalar *v) 79b49de8d1SLois Curfman McInnes { 80b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 81b49de8d1SLois Curfman McInnes int ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row; 82b49de8d1SLois Curfman McInnes 833a40ed3dSBarry Smith PetscFunctionBegin; 84b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 8529bbc08cSBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 86273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 87b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 88b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 89b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 9029bbc08cSBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 91273d9f13SBarry Smith if (idxn[j] >= mat->N) { 9229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 93a8c6a408SBarry Smith } 94b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 95b49de8d1SLois Curfman McInnes } 96a8c6a408SBarry Smith } else { 9729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 988965ea79SLois Curfman McInnes } 998965ea79SLois Curfman McInnes } 1003a40ed3dSBarry Smith PetscFunctionReturn(0); 1018965ea79SLois Curfman McInnes } 1028965ea79SLois Curfman McInnes 1034a2ae208SSatish Balay #undef __FUNCT__ 1044a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 105*87828ca2SBarry Smith int MatGetArray_MPIDense(Mat A,PetscScalar **array) 106ff14e315SSatish Balay { 107ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 108ff14e315SSatish Balay int ierr; 109ff14e315SSatish Balay 1103a40ed3dSBarry Smith PetscFunctionBegin; 111ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1123a40ed3dSBarry Smith PetscFunctionReturn(0); 113ff14e315SSatish Balay } 114ff14e315SSatish Balay 1154a2ae208SSatish Balay #undef __FUNCT__ 1164a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 117ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 118ca3fa75bSLois Curfman McInnes { 119ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 120ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 121cfce73b9SSatish Balay int i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 122*87828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 123ca3fa75bSLois Curfman McInnes Mat newmat; 124ca3fa75bSLois Curfman McInnes 125ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 126ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 127ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 128b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 129b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 130ca3fa75bSLois Curfman McInnes 131ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1327eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1337eba5e9cSLois Curfman McInnes original matrix! */ 134ca3fa75bSLois Curfman McInnes 135ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1367eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 137ca3fa75bSLois Curfman McInnes 138ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 139ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 14029bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 1417eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 142ca3fa75bSLois Curfman McInnes newmat = *B; 143ca3fa75bSLois Curfman McInnes } else { 144ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 14532828cfdSBarry Smith ierr = MatCreateMPIDense(A->comm,nrows,cs,PETSC_DECIDE,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 146ca3fa75bSLois Curfman McInnes } 147ca3fa75bSLois Curfman McInnes 148ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 149ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 150ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 151ca3fa75bSLois Curfman McInnes 152ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 153ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 154ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 1557eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 156ca3fa75bSLois Curfman McInnes } 157ca3fa75bSLois Curfman McInnes } 158ca3fa75bSLois Curfman McInnes 159ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 160ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 161ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 162ca3fa75bSLois Curfman McInnes 163ca3fa75bSLois Curfman McInnes /* Free work space */ 164ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 165ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 166ca3fa75bSLois Curfman McInnes *B = newmat; 167ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 168ca3fa75bSLois Curfman McInnes } 169ca3fa75bSLois Curfman McInnes 1704a2ae208SSatish Balay #undef __FUNCT__ 1714a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 172*87828ca2SBarry Smith int MatRestoreArray_MPIDense(Mat A,PetscScalar **array) 173ff14e315SSatish Balay { 1743a40ed3dSBarry Smith PetscFunctionBegin; 1753a40ed3dSBarry Smith PetscFunctionReturn(0); 176ff14e315SSatish Balay } 177ff14e315SSatish Balay 1784a2ae208SSatish Balay #undef __FUNCT__ 1794a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 1808f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 1818965ea79SLois Curfman McInnes { 18239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 1838965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 184d36fbae8SSatish Balay int ierr,nstash,reallocs; 1858965ea79SLois Curfman McInnes InsertMode addv; 1868965ea79SLois Curfman McInnes 1873a40ed3dSBarry Smith PetscFunctionBegin; 1888965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 189ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1907056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 19129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 1928965ea79SLois Curfman McInnes } 193e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1948965ea79SLois Curfman McInnes 1958798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 1968798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 197b0a32e0cSBarry Smith PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs); 1983a40ed3dSBarry Smith PetscFunctionReturn(0); 1998965ea79SLois Curfman McInnes } 2008965ea79SLois Curfman McInnes 2014a2ae208SSatish Balay #undef __FUNCT__ 2024a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 2038f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2048965ea79SLois Curfman McInnes { 20539ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2067ef1d9bdSSatish Balay int i,n,ierr,*row,*col,flg,j,rstart,ncols; 207*87828ca2SBarry Smith PetscScalar *val; 208e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2098965ea79SLois Curfman McInnes 2103a40ed3dSBarry Smith PetscFunctionBegin; 2118965ea79SLois Curfman McInnes /* wait on receives */ 2127ef1d9bdSSatish Balay while (1) { 2138798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2147ef1d9bdSSatish Balay if (!flg) break; 2158965ea79SLois Curfman McInnes 2167ef1d9bdSSatish Balay for (i=0; i<n;) { 2177ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2187ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2197ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2207ef1d9bdSSatish Balay else ncols = n-i; 2217ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2227ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2237ef1d9bdSSatish Balay i = j; 2248965ea79SLois Curfman McInnes } 2257ef1d9bdSSatish Balay } 2268798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2278965ea79SLois Curfman McInnes 22839ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 22939ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2308965ea79SLois Curfman McInnes 2316d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 23239ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2338965ea79SLois Curfman McInnes } 2343a40ed3dSBarry Smith PetscFunctionReturn(0); 2358965ea79SLois Curfman McInnes } 2368965ea79SLois Curfman McInnes 2374a2ae208SSatish Balay #undef __FUNCT__ 2384a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 2398f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A) 2408965ea79SLois Curfman McInnes { 2413a40ed3dSBarry Smith int ierr; 24239ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2433a40ed3dSBarry Smith 2443a40ed3dSBarry Smith PetscFunctionBegin; 2453a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 2478965ea79SLois Curfman McInnes } 2488965ea79SLois Curfman McInnes 2494a2ae208SSatish Balay #undef __FUNCT__ 2504a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIDense" 2518f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs) 2524e220ebcSLois Curfman McInnes { 2533a40ed3dSBarry Smith PetscFunctionBegin; 2544e220ebcSLois Curfman McInnes *bs = 1; 2553a40ed3dSBarry Smith PetscFunctionReturn(0); 2564e220ebcSLois Curfman McInnes } 2574e220ebcSLois Curfman McInnes 2588965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2598965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2608965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2613501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2628965ea79SLois Curfman McInnes routine. 2638965ea79SLois Curfman McInnes */ 2644a2ae208SSatish Balay #undef __FUNCT__ 2654a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 266*87828ca2SBarry Smith int MatZeroRows_MPIDense(Mat A,IS is,PetscScalar *diag) 2678965ea79SLois Curfman McInnes { 26839ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2698965ea79SLois Curfman McInnes int i,ierr,N,*rows,*owners = l->rowners,size = l->size; 27035d8aa7fSBarry Smith int *procs,*nprocs,j,idx,nsends,*work; 2718965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2728965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2738965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2748965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2758965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2768965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2778965ea79SLois Curfman McInnes IS istmp; 27835d8aa7fSBarry Smith PetscTruth found; 2798965ea79SLois Curfman McInnes 2803a40ed3dSBarry Smith PetscFunctionBegin; 281b9b97703SBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 2828965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 2838965ea79SLois Curfman McInnes 2848965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 285b0a32e0cSBarry Smith ierr = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr); 286549d3d68SSatish Balay ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 287549d3d68SSatish Balay procs = nprocs + size; 288b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/ 2898965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 2908965ea79SLois Curfman McInnes idx = rows[i]; 29135d8aa7fSBarry Smith found = PETSC_FALSE; 2928965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 2938965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 29435d8aa7fSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = PETSC_TRUE; break; 2958965ea79SLois Curfman McInnes } 2968965ea79SLois Curfman McInnes } 29729bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 2988965ea79SLois Curfman McInnes } 2998965ea79SLois Curfman McInnes nsends = 0; for (i=0; i<size; i++) { nsends += procs[i];} 3008965ea79SLois Curfman McInnes 3018965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 302b0a32e0cSBarry Smith ierr = PetscMalloc(2*size*sizeof(int),&work);CHKERRQ(ierr); 3036831982aSBarry Smith ierr = MPI_Allreduce(nprocs,work,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr); 3048965ea79SLois Curfman McInnes nmax = work[rank]; 3056831982aSBarry Smith nrecvs = work[size+rank]; 306606d414cSSatish Balay ierr = PetscFree(work);CHKERRQ(ierr); 3078965ea79SLois Curfman McInnes 3088965ea79SLois Curfman McInnes /* post receives: */ 309b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr); 310b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3118965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 312ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3138965ea79SLois Curfman McInnes } 3148965ea79SLois Curfman McInnes 3158965ea79SLois Curfman McInnes /* do sends: 3168965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3178965ea79SLois Curfman McInnes the ith processor 3188965ea79SLois Curfman McInnes */ 319b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr); 320b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 321b0a32e0cSBarry Smith ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr); 3228965ea79SLois Curfman McInnes starts[0] = 0; 3238965ea79SLois Curfman McInnes for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];} 3248965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3258965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3268965ea79SLois Curfman McInnes } 3278965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3288965ea79SLois Curfman McInnes 3298965ea79SLois Curfman McInnes starts[0] = 0; 3308965ea79SLois Curfman McInnes for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[i-1];} 3318965ea79SLois Curfman McInnes count = 0; 3328965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 3338965ea79SLois Curfman McInnes if (procs[i]) { 334ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3358965ea79SLois Curfman McInnes } 3368965ea79SLois Curfman McInnes } 337606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3388965ea79SLois Curfman McInnes 3398965ea79SLois Curfman McInnes base = owners[rank]; 3408965ea79SLois Curfman McInnes 3418965ea79SLois Curfman McInnes /* wait on receives */ 342b0a32e0cSBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr); 3438965ea79SLois Curfman McInnes source = lens + nrecvs; 3448965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3458965ea79SLois Curfman McInnes while (count) { 346ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3478965ea79SLois Curfman McInnes /* unpack receives into our local space */ 348ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 3498965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3508965ea79SLois Curfman McInnes lens[imdex] = n; 3518965ea79SLois Curfman McInnes slen += n; 3528965ea79SLois Curfman McInnes count--; 3538965ea79SLois Curfman McInnes } 354606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3558965ea79SLois Curfman McInnes 3568965ea79SLois Curfman McInnes /* move the data into the send scatter */ 357b0a32e0cSBarry Smith ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr); 3588965ea79SLois Curfman McInnes count = 0; 3598965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 3608965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3618965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 3628965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3638965ea79SLois Curfman McInnes } 3648965ea79SLois Curfman McInnes } 365606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 366606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 367606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 368606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3698965ea79SLois Curfman McInnes 3708965ea79SLois Curfman McInnes /* actually zap the local rows */ 371029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 372b0a32e0cSBarry Smith PetscLogObjectParent(A,istmp); 373606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3748965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 3758965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 3768965ea79SLois Curfman McInnes 3778965ea79SLois Curfman McInnes /* wait on sends */ 3788965ea79SLois Curfman McInnes if (nsends) { 379b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 380ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 381606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 3828965ea79SLois Curfman McInnes } 383606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 384606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 3858965ea79SLois Curfman McInnes 3863a40ed3dSBarry Smith PetscFunctionReturn(0); 3878965ea79SLois Curfman McInnes } 3888965ea79SLois Curfman McInnes 3894a2ae208SSatish Balay #undef __FUNCT__ 3904a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 3918f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3928965ea79SLois Curfman McInnes { 39339ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 3948965ea79SLois Curfman McInnes int ierr; 395c456f294SBarry Smith 3963a40ed3dSBarry Smith PetscFunctionBegin; 39743a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39843a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39944cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4003a40ed3dSBarry Smith PetscFunctionReturn(0); 4018965ea79SLois Curfman McInnes } 4028965ea79SLois Curfman McInnes 4034a2ae208SSatish Balay #undef __FUNCT__ 4044a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 4058f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4068965ea79SLois Curfman McInnes { 40739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4088965ea79SLois Curfman McInnes int ierr; 409c456f294SBarry Smith 4103a40ed3dSBarry Smith PetscFunctionBegin; 41143a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41243a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41344cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4143a40ed3dSBarry Smith PetscFunctionReturn(0); 4158965ea79SLois Curfman McInnes } 4168965ea79SLois Curfman McInnes 4174a2ae208SSatish Balay #undef __FUNCT__ 4184a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 4197c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 420096963f5SLois Curfman McInnes { 421096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 422096963f5SLois Curfman McInnes int ierr; 423*87828ca2SBarry Smith PetscScalar zero = 0.0; 424096963f5SLois Curfman McInnes 4253a40ed3dSBarry Smith PetscFunctionBegin; 4263501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 4277c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 428537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 429537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4303a40ed3dSBarry Smith PetscFunctionReturn(0); 431096963f5SLois Curfman McInnes } 432096963f5SLois Curfman McInnes 4334a2ae208SSatish Balay #undef __FUNCT__ 4344a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 4357c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 436096963f5SLois Curfman McInnes { 437096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 438096963f5SLois Curfman McInnes int ierr; 439096963f5SLois Curfman McInnes 4403a40ed3dSBarry Smith PetscFunctionBegin; 4413501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4427c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 443537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 444537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4453a40ed3dSBarry Smith PetscFunctionReturn(0); 446096963f5SLois Curfman McInnes } 447096963f5SLois Curfman McInnes 4484a2ae208SSatish Balay #undef __FUNCT__ 4494a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 4508f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 4518965ea79SLois Curfman McInnes { 45239ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 453096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 454273d9f13SBarry Smith int ierr,len,i,n,m = A->m,radd; 455*87828ca2SBarry Smith PetscScalar *x,zero = 0.0; 456ed3cc1f0SBarry Smith 4573a40ed3dSBarry Smith PetscFunctionBegin; 458273d9f13SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 459096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x);CHKERRQ(ierr); 460096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 461273d9f13SBarry Smith if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 462273d9f13SBarry Smith len = PetscMin(a->A->m,a->A->n); 4637ddc982cSLois Curfman McInnes radd = a->rstart*m; 46444cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 465096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 466096963f5SLois Curfman McInnes } 4679a8c540fSSatish Balay ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 4683a40ed3dSBarry Smith PetscFunctionReturn(0); 4698965ea79SLois Curfman McInnes } 4708965ea79SLois Curfman McInnes 4714a2ae208SSatish Balay #undef __FUNCT__ 4724a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 473e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 4748965ea79SLois Curfman McInnes { 4753501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4768965ea79SLois Curfman McInnes int ierr; 477ed3cc1f0SBarry Smith 4783a40ed3dSBarry Smith PetscFunctionBegin; 47994d884c6SBarry Smith 480aa482453SBarry Smith #if defined(PETSC_USE_LOG) 481b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N); 4828965ea79SLois Curfman McInnes #endif 4838798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 484606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 4853501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 4863501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4873501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 488622d7880SLois Curfman McInnes if (mdn->factor) { 489606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 490606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 491606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 492606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 493622d7880SLois Curfman McInnes } 494606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 4953a40ed3dSBarry Smith PetscFunctionReturn(0); 4968965ea79SLois Curfman McInnes } 49739ddd567SLois Curfman McInnes 4984a2ae208SSatish Balay #undef __FUNCT__ 4994a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 500b0a32e0cSBarry Smith static int MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5018965ea79SLois Curfman McInnes { 50239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 5038965ea79SLois Curfman McInnes int ierr; 5047056b6fcSBarry Smith 5053a40ed3dSBarry Smith PetscFunctionBegin; 50639ddd567SLois Curfman McInnes if (mdn->size == 1) { 50739ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5088965ea79SLois Curfman McInnes } 50929bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5103a40ed3dSBarry Smith PetscFunctionReturn(0); 5118965ea79SLois Curfman McInnes } 5128965ea79SLois Curfman McInnes 5134a2ae208SSatish Balay #undef __FUNCT__ 5144a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 515b0a32e0cSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5168965ea79SLois Curfman McInnes { 51739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 518fb9695e5SSatish Balay int ierr,size = mdn->size,rank = mdn->rank; 519b0a32e0cSBarry Smith PetscViewerType vtype; 520f1af5d2fSBarry Smith PetscTruth isascii,isdraw; 521b0a32e0cSBarry Smith PetscViewer sviewer; 522f3ef73ceSBarry Smith PetscViewerFormat format; 5238965ea79SLois Curfman McInnes 5243a40ed3dSBarry Smith PetscFunctionBegin; 525b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 526fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 527f1af5d2fSBarry Smith if (isascii) { 528b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 529b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 530fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO_LONG) { 5314e220ebcSLois Curfman McInnes MatInfo info; 532888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 533b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m, 5346831982aSBarry Smith (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr); 535b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5363501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5373a40ed3dSBarry Smith PetscFunctionReturn(0); 538fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 5393a40ed3dSBarry Smith PetscFunctionReturn(0); 5408965ea79SLois Curfman McInnes } 541f1af5d2fSBarry Smith } else if (isdraw) { 542b0a32e0cSBarry Smith PetscDraw draw; 543f1af5d2fSBarry Smith PetscTruth isnull; 544f1af5d2fSBarry Smith 545b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 546b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 547f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 548f1af5d2fSBarry Smith } 54977ed5343SBarry Smith 5508965ea79SLois Curfman McInnes if (size == 1) { 55139ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5523a40ed3dSBarry Smith } else { 5538965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5548965ea79SLois Curfman McInnes Mat A; 555273d9f13SBarry Smith int M = mat->M,N = mat->N,m,row,i,nz,*cols; 556*87828ca2SBarry Smith PetscScalar *vals; 5578965ea79SLois Curfman McInnes 5588965ea79SLois Curfman McInnes if (!rank) { 559f1af5d2fSBarry Smith ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 5603a40ed3dSBarry Smith } else { 561f1af5d2fSBarry Smith ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 5628965ea79SLois Curfman McInnes } 563b0a32e0cSBarry Smith PetscLogObjectParent(mat,A); 5648965ea79SLois Curfman McInnes 56539ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 56639ddd567SLois Curfman McInnes but it's quick for now */ 567273d9f13SBarry Smith row = mdn->rstart; m = mdn->A->m; 5688965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 56939ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 57039ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 57139ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 57239ddd567SLois Curfman McInnes row++; 5738965ea79SLois Curfman McInnes } 5748965ea79SLois Curfman McInnes 5756d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5766d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 577b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 578b9b97703SBarry Smith if (!rank) { 5796831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 5808965ea79SLois Curfman McInnes } 581b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 582b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5838965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 5848965ea79SLois Curfman McInnes } 5853a40ed3dSBarry Smith PetscFunctionReturn(0); 5868965ea79SLois Curfman McInnes } 5878965ea79SLois Curfman McInnes 5884a2ae208SSatish Balay #undef __FUNCT__ 5894a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 590b0a32e0cSBarry Smith int MatView_MPIDense(Mat mat,PetscViewer viewer) 5918965ea79SLois Curfman McInnes { 59239ddd567SLois Curfman McInnes int ierr; 593f1af5d2fSBarry Smith PetscTruth isascii,isbinary,isdraw,issocket; 5948965ea79SLois Curfman McInnes 595433994e6SBarry Smith PetscFunctionBegin; 5960f5bd95cSBarry Smith 597b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 598fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 599b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 600fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 6010f5bd95cSBarry Smith 602f1af5d2fSBarry Smith if (isascii || issocket || isdraw) { 603f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 6040f5bd95cSBarry Smith } else if (isbinary) { 6053a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6065cd90555SBarry Smith } else { 60729bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6088965ea79SLois Curfman McInnes } 6093a40ed3dSBarry Smith PetscFunctionReturn(0); 6108965ea79SLois Curfman McInnes } 6118965ea79SLois Curfman McInnes 6124a2ae208SSatish Balay #undef __FUNCT__ 6134a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 6148f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6158965ea79SLois Curfman McInnes { 6163501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 6173501a2bdSLois Curfman McInnes Mat mdn = mat->A; 6184e220ebcSLois Curfman McInnes int ierr; 619329f5518SBarry Smith PetscReal isend[5],irecv[5]; 6208965ea79SLois Curfman McInnes 6213a40ed3dSBarry Smith PetscFunctionBegin; 622273d9f13SBarry Smith info->rows_global = (double)A->M; 623273d9f13SBarry Smith info->columns_global = (double)A->N; 624273d9f13SBarry Smith info->rows_local = (double)A->m; 625273d9f13SBarry Smith info->columns_local = (double)A->N; 6264e220ebcSLois Curfman McInnes info->block_size = 1.0; 6274e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6284e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6294e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6308965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6314e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6324e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6334e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6344e220ebcSLois Curfman McInnes info->memory = isend[3]; 6354e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6368965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 637f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 6384e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6394e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6404e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6414e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6424e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6438965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 644f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 6454e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6464e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6474e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6484e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6494e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6508965ea79SLois Curfman McInnes } 6514e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6524e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6534e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6543a40ed3dSBarry Smith PetscFunctionReturn(0); 6558965ea79SLois Curfman McInnes } 6568965ea79SLois Curfman McInnes 6574a2ae208SSatish Balay #undef __FUNCT__ 6584a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 6598f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 6608965ea79SLois Curfman McInnes { 66139ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 662273d9f13SBarry Smith int ierr; 6638965ea79SLois Curfman McInnes 6643a40ed3dSBarry Smith PetscFunctionBegin; 66512c028f9SKris Buschelman switch (op) { 66612c028f9SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 66712c028f9SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 66812c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 66912c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 67012c028f9SKris Buschelman case MAT_COLUMNS_SORTED: 67112c028f9SKris Buschelman case MAT_COLUMNS_UNSORTED: 672273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 67312c028f9SKris Buschelman break; 67412c028f9SKris Buschelman case MAT_ROW_ORIENTED: 675273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 676273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 67712c028f9SKris Buschelman break; 67812c028f9SKris Buschelman case MAT_ROWS_SORTED: 67912c028f9SKris Buschelman case MAT_ROWS_UNSORTED: 68012c028f9SKris Buschelman case MAT_YES_NEW_DIAGONALS: 68112c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 682d03495bdSKris Buschelman case MAT_USE_SINGLE_PRECISION_SOLVES: 683b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 68412c028f9SKris Buschelman break; 68512c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 686273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 687273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 68812c028f9SKris Buschelman break; 68912c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 690273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 69112c028f9SKris Buschelman break; 69212c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 69329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 69412c028f9SKris Buschelman break; 69512c028f9SKris Buschelman default: 69629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 69712c028f9SKris Buschelman break; 6983a40ed3dSBarry Smith } 6993a40ed3dSBarry Smith PetscFunctionReturn(0); 7008965ea79SLois Curfman McInnes } 7018965ea79SLois Curfman McInnes 7024a2ae208SSatish Balay #undef __FUNCT__ 7034a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_MPIDense" 7048f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 7058965ea79SLois Curfman McInnes { 7063501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7073a40ed3dSBarry Smith 7083a40ed3dSBarry Smith PetscFunctionBegin; 7094c49b128SBarry Smith if (m) *m = mat->rstart; 7104c49b128SBarry Smith if (n) *n = mat->rend; 7113a40ed3dSBarry Smith PetscFunctionReturn(0); 7128965ea79SLois Curfman McInnes } 7138965ea79SLois Curfman McInnes 7144a2ae208SSatish Balay #undef __FUNCT__ 7154a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIDense" 716*87828ca2SBarry Smith int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v) 7178965ea79SLois Curfman McInnes { 7183501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7193a40ed3dSBarry Smith int lrow,rstart = mat->rstart,rend = mat->rend,ierr; 7208965ea79SLois Curfman McInnes 7213a40ed3dSBarry Smith PetscFunctionBegin; 72229bbc08cSBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 7238965ea79SLois Curfman McInnes lrow = row - rstart; 7243a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 7253a40ed3dSBarry Smith PetscFunctionReturn(0); 7268965ea79SLois Curfman McInnes } 7278965ea79SLois Curfman McInnes 7284a2ae208SSatish Balay #undef __FUNCT__ 7294a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIDense" 730*87828ca2SBarry Smith int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v) 7318965ea79SLois Curfman McInnes { 732606d414cSSatish Balay int ierr; 733606d414cSSatish Balay 7343a40ed3dSBarry Smith PetscFunctionBegin; 735606d414cSSatish Balay if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 736606d414cSSatish Balay if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 7373a40ed3dSBarry Smith PetscFunctionReturn(0); 7388965ea79SLois Curfman McInnes } 7398965ea79SLois Curfman McInnes 7404a2ae208SSatish Balay #undef __FUNCT__ 7414a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 7425b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7435b2fa520SLois Curfman McInnes { 7445b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7455b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 746*87828ca2SBarry Smith PetscScalar *l,*r,x,*v; 747273d9f13SBarry Smith int ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 7485b2fa520SLois Curfman McInnes 7495b2fa520SLois Curfman McInnes PetscFunctionBegin; 75072d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7515b2fa520SLois Curfman McInnes if (ll) { 75272d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 75329bbc08cSBarry Smith if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 7545b2fa520SLois Curfman McInnes ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7555b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7565b2fa520SLois Curfman McInnes x = l[i]; 7575b2fa520SLois Curfman McInnes v = mat->v + i; 7585b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7595b2fa520SLois Curfman McInnes } 7605b2fa520SLois Curfman McInnes ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 761b0a32e0cSBarry Smith PetscLogFlops(n*m); 7625b2fa520SLois Curfman McInnes } 7635b2fa520SLois Curfman McInnes if (rr) { 76472d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 76529bbc08cSBarry Smith if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 7665b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7675b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7685b2fa520SLois Curfman McInnes ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 7695b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 7705b2fa520SLois Curfman McInnes x = r[i]; 7715b2fa520SLois Curfman McInnes v = mat->v + i*m; 7725b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 7735b2fa520SLois Curfman McInnes } 77472d926a5SLois Curfman McInnes ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 775b0a32e0cSBarry Smith PetscLogFlops(n*m); 7765b2fa520SLois Curfman McInnes } 7775b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 7785b2fa520SLois Curfman McInnes } 7795b2fa520SLois Curfman McInnes 7804a2ae208SSatish Balay #undef __FUNCT__ 7814a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 782329f5518SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *norm) 783096963f5SLois Curfman McInnes { 7843501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7853501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 7863501a2bdSLois Curfman McInnes int ierr,i,j; 787329f5518SBarry Smith PetscReal sum = 0.0; 788*87828ca2SBarry Smith PetscScalar *v = mat->v; 7893501a2bdSLois Curfman McInnes 7903a40ed3dSBarry Smith PetscFunctionBegin; 7913501a2bdSLois Curfman McInnes if (mdn->size == 1) { 7923501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm);CHKERRQ(ierr); 7933501a2bdSLois Curfman McInnes } else { 7943501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 795273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 796aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 797329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 7983501a2bdSLois Curfman McInnes #else 7993501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 8003501a2bdSLois Curfman McInnes #endif 8013501a2bdSLois Curfman McInnes } 802ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 8033501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 804b0a32e0cSBarry Smith PetscLogFlops(2*mdn->A->n*mdn->A->m); 8053a40ed3dSBarry Smith } else if (type == NORM_1) { 806329f5518SBarry Smith PetscReal *tmp,*tmp2; 807b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 808273d9f13SBarry Smith tmp2 = tmp + A->N; 809273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 810096963f5SLois Curfman McInnes *norm = 0.0; 8113501a2bdSLois Curfman McInnes v = mat->v; 812273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 813273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 81467e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8153501a2bdSLois Curfman McInnes } 8163501a2bdSLois Curfman McInnes } 817273d9f13SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 818273d9f13SBarry Smith for (j=0; j<A->N; j++) { 8193501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 8203501a2bdSLois Curfman McInnes } 821606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 822b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 8233a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 824329f5518SBarry Smith PetscReal ntemp; 8253501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 826ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 8273a40ed3dSBarry Smith } else { 82829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8293501a2bdSLois Curfman McInnes } 8303501a2bdSLois Curfman McInnes } 8313a40ed3dSBarry Smith PetscFunctionReturn(0); 8323501a2bdSLois Curfman McInnes } 8333501a2bdSLois Curfman McInnes 8344a2ae208SSatish Balay #undef __FUNCT__ 8354a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 8368f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 8373501a2bdSLois Curfman McInnes { 8383501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8393501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8403501a2bdSLois Curfman McInnes Mat B; 841273d9f13SBarry Smith int M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 8423501a2bdSLois Curfman McInnes int j,i,ierr; 843*87828ca2SBarry Smith PetscScalar *v; 8443501a2bdSLois Curfman McInnes 8453a40ed3dSBarry Smith PetscFunctionBegin; 8467c922b88SBarry Smith if (!matout && M != N) { 84729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 8487056b6fcSBarry Smith } 8497056b6fcSBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 8503501a2bdSLois Curfman McInnes 851273d9f13SBarry Smith m = a->A->m; n = a->A->n; v = Aloc->v; 852b0a32e0cSBarry Smith ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr); 8533501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 8543501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8553501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8563501a2bdSLois Curfman McInnes v += m; 8573501a2bdSLois Curfman McInnes } 858606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 8596d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8606d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8617c922b88SBarry Smith if (matout) { 8623501a2bdSLois Curfman McInnes *matout = B; 8633501a2bdSLois Curfman McInnes } else { 864273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 8653501a2bdSLois Curfman McInnes } 8663a40ed3dSBarry Smith PetscFunctionReturn(0); 867096963f5SLois Curfman McInnes } 868096963f5SLois Curfman McInnes 869d9eff348SSatish Balay #include "petscblaslapack.h" 8704a2ae208SSatish Balay #undef __FUNCT__ 8714a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 872*87828ca2SBarry Smith int MatScale_MPIDense(PetscScalar *alpha,Mat inA) 87344cd7ae7SLois Curfman McInnes { 87444cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 87544cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 87644cd7ae7SLois Curfman McInnes int one = 1,nz; 87744cd7ae7SLois Curfman McInnes 8783a40ed3dSBarry Smith PetscFunctionBegin; 879273d9f13SBarry Smith nz = inA->m*inA->N; 88044cd7ae7SLois Curfman McInnes BLscal_(&nz,alpha,a->v,&one); 881b0a32e0cSBarry Smith PetscLogFlops(nz); 8823a40ed3dSBarry Smith PetscFunctionReturn(0); 88344cd7ae7SLois Curfman McInnes } 88444cd7ae7SLois Curfman McInnes 8855609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 886ca44d042SBarry Smith EXTERN int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **); 8878965ea79SLois Curfman McInnes 8884a2ae208SSatish Balay #undef __FUNCT__ 8894a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 890273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A) 891273d9f13SBarry Smith { 892273d9f13SBarry Smith int ierr; 893273d9f13SBarry Smith 894273d9f13SBarry Smith PetscFunctionBegin; 895273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 896273d9f13SBarry Smith PetscFunctionReturn(0); 897273d9f13SBarry Smith } 898273d9f13SBarry Smith 8998965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 90009dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 90109dc0095SBarry Smith MatGetRow_MPIDense, 90209dc0095SBarry Smith MatRestoreRow_MPIDense, 90309dc0095SBarry Smith MatMult_MPIDense, 90409dc0095SBarry Smith MatMultAdd_MPIDense, 9057c922b88SBarry Smith MatMultTranspose_MPIDense, 9067c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 9078965ea79SLois Curfman McInnes 0, 90809dc0095SBarry Smith 0, 90909dc0095SBarry Smith 0, 91009dc0095SBarry Smith 0, 91109dc0095SBarry Smith 0, 91209dc0095SBarry Smith 0, 91309dc0095SBarry Smith 0, 91409dc0095SBarry Smith MatTranspose_MPIDense, 91509dc0095SBarry Smith MatGetInfo_MPIDense,0, 91609dc0095SBarry Smith MatGetDiagonal_MPIDense, 9175b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 91809dc0095SBarry Smith MatNorm_MPIDense, 91909dc0095SBarry Smith MatAssemblyBegin_MPIDense, 92009dc0095SBarry Smith MatAssemblyEnd_MPIDense, 92109dc0095SBarry Smith 0, 92209dc0095SBarry Smith MatSetOption_MPIDense, 92309dc0095SBarry Smith MatZeroEntries_MPIDense, 92409dc0095SBarry Smith MatZeroRows_MPIDense, 92509dc0095SBarry Smith 0, 92609dc0095SBarry Smith 0, 92709dc0095SBarry Smith 0, 92809dc0095SBarry Smith 0, 929273d9f13SBarry Smith MatSetUpPreallocation_MPIDense, 930273d9f13SBarry Smith 0, 93139ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 93209dc0095SBarry Smith 0, 93309dc0095SBarry Smith 0, 93409dc0095SBarry Smith MatGetArray_MPIDense, 93509dc0095SBarry Smith MatRestoreArray_MPIDense, 9365609ef8eSBarry Smith MatDuplicate_MPIDense, 93709dc0095SBarry Smith 0, 93809dc0095SBarry Smith 0, 93909dc0095SBarry Smith 0, 94009dc0095SBarry Smith 0, 94109dc0095SBarry Smith 0, 9422ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 94309dc0095SBarry Smith 0, 94409dc0095SBarry Smith MatGetValues_MPIDense, 94509dc0095SBarry Smith 0, 94609dc0095SBarry Smith 0, 94709dc0095SBarry Smith MatScale_MPIDense, 94809dc0095SBarry Smith 0, 94909dc0095SBarry Smith 0, 95009dc0095SBarry Smith 0, 95109dc0095SBarry Smith MatGetBlockSize_MPIDense, 95209dc0095SBarry Smith 0, 95309dc0095SBarry Smith 0, 95409dc0095SBarry Smith 0, 95509dc0095SBarry Smith 0, 95609dc0095SBarry Smith 0, 95709dc0095SBarry Smith 0, 95809dc0095SBarry Smith 0, 95909dc0095SBarry Smith 0, 96009dc0095SBarry Smith 0, 961ca3fa75bSLois Curfman McInnes MatGetSubMatrix_MPIDense, 962b9b97703SBarry Smith MatDestroy_MPIDense, 963b9b97703SBarry Smith MatView_MPIDense, 9648a124369SBarry Smith MatGetPetscMaps_Petsc}; 9658965ea79SLois Curfman McInnes 966273d9f13SBarry Smith EXTERN_C_BEGIN 9674a2ae208SSatish Balay #undef __FUNCT__ 9684a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 969273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat) 970273d9f13SBarry Smith { 971273d9f13SBarry Smith Mat_MPIDense *a; 972273d9f13SBarry Smith int ierr,i; 973273d9f13SBarry Smith 974273d9f13SBarry Smith PetscFunctionBegin; 975b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 976b0a32e0cSBarry Smith mat->data = (void*)a; 977273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 978273d9f13SBarry Smith mat->factor = 0; 979273d9f13SBarry Smith mat->mapping = 0; 980273d9f13SBarry Smith 981273d9f13SBarry Smith a->factor = 0; 982273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 983273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 984273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 985273d9f13SBarry Smith 986273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 987273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 988273d9f13SBarry Smith a->nvec = mat->n; 989273d9f13SBarry Smith 990273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 991273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 9928a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 9938a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 994273d9f13SBarry Smith 995273d9f13SBarry Smith /* build local table of row and column ownerships */ 996b0a32e0cSBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr); 997273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 998b0a32e0cSBarry Smith PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 999273d9f13SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1000273d9f13SBarry Smith a->rowners[0] = 0; 1001273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1002273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1003273d9f13SBarry Smith } 1004273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1005273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 1006273d9f13SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1007273d9f13SBarry Smith a->cowners[0] = 0; 1008273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1009273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1010273d9f13SBarry Smith } 1011273d9f13SBarry Smith 1012273d9f13SBarry Smith /* build cache for off array entries formed */ 1013273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1014273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1015273d9f13SBarry Smith 1016273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1017273d9f13SBarry Smith a->lvec = 0; 1018273d9f13SBarry Smith a->Mvctx = 0; 1019273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1020273d9f13SBarry Smith 1021273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1022273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1023273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1024273d9f13SBarry Smith PetscFunctionReturn(0); 1025273d9f13SBarry Smith } 1026273d9f13SBarry Smith EXTERN_C_END 1027273d9f13SBarry Smith 10284a2ae208SSatish Balay #undef __FUNCT__ 10294a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1030273d9f13SBarry Smith /*@C 1031273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1032273d9f13SBarry Smith 1033273d9f13SBarry Smith Not collective 1034273d9f13SBarry Smith 1035273d9f13SBarry Smith Input Parameters: 1036273d9f13SBarry Smith . A - the matrix 1037273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1038273d9f13SBarry Smith to control all matrix memory allocation. 1039273d9f13SBarry Smith 1040273d9f13SBarry Smith Notes: 1041273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1042273d9f13SBarry Smith storage by columns. 1043273d9f13SBarry Smith 1044273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1045273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1046273d9f13SBarry Smith set data=PETSC_NULL. 1047273d9f13SBarry Smith 1048273d9f13SBarry Smith Level: intermediate 1049273d9f13SBarry Smith 1050273d9f13SBarry Smith .keywords: matrix,dense, parallel 1051273d9f13SBarry Smith 1052273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1053273d9f13SBarry Smith @*/ 1054*87828ca2SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1055273d9f13SBarry Smith { 1056273d9f13SBarry Smith Mat_MPIDense *a; 1057273d9f13SBarry Smith int ierr; 1058273d9f13SBarry Smith PetscTruth flg2; 1059273d9f13SBarry Smith 1060273d9f13SBarry Smith PetscFunctionBegin; 1061273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATMPIDENSE,&flg2);CHKERRQ(ierr); 1062273d9f13SBarry Smith if (!flg2) PetscFunctionReturn(0); 1063273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 1064273d9f13SBarry Smith /* Note: For now, when data is specified above, this assumes the user correctly 1065273d9f13SBarry Smith allocates the local dense storage space. We should add error checking. */ 1066273d9f13SBarry Smith 1067273d9f13SBarry Smith a = (Mat_MPIDense*)mat->data; 1068273d9f13SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,mat->m,mat->N,data,&a->A);CHKERRQ(ierr); 1069b0a32e0cSBarry Smith PetscLogObjectParent(mat,a->A); 1070273d9f13SBarry Smith PetscFunctionReturn(0); 1071273d9f13SBarry Smith } 1072273d9f13SBarry Smith 10734a2ae208SSatish Balay #undef __FUNCT__ 10744a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 10758965ea79SLois Curfman McInnes /*@C 107639ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 10778965ea79SLois Curfman McInnes 1078db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1079db81eaa0SLois Curfman McInnes 10808965ea79SLois Curfman McInnes Input Parameters: 1081db81eaa0SLois Curfman McInnes + comm - MPI communicator 10828965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1083db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 10848965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1085db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 1086db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1087dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 10888965ea79SLois Curfman McInnes 10898965ea79SLois Curfman McInnes Output Parameter: 1090477f1c0bSLois Curfman McInnes . A - the matrix 10918965ea79SLois Curfman McInnes 1092b259b22eSLois Curfman McInnes Notes: 109339ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 109439ddd567SLois Curfman McInnes storage by columns. 10958965ea79SLois Curfman McInnes 109618f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 109718f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1098b4fd4287SBarry Smith set data=PETSC_NULL. 109918f449edSLois Curfman McInnes 11008965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 11018965ea79SLois Curfman McInnes (possibly both). 11028965ea79SLois Curfman McInnes 1103027ccd11SLois Curfman McInnes Level: intermediate 1104027ccd11SLois Curfman McInnes 110539ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 11068965ea79SLois Curfman McInnes 110739ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 11088965ea79SLois Curfman McInnes @*/ 1109*87828ca2SBarry Smith int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A) 11108965ea79SLois Curfman McInnes { 1111273d9f13SBarry Smith int ierr,size; 11128965ea79SLois Curfman McInnes 11133a40ed3dSBarry Smith PetscFunctionBegin; 1114273d9f13SBarry Smith ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr); 1115273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1116273d9f13SBarry Smith if (size > 1) { 1117273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1118273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1119273d9f13SBarry Smith } else { 1120273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1121273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 11228c469469SLois Curfman McInnes } 11233a40ed3dSBarry Smith PetscFunctionReturn(0); 11248965ea79SLois Curfman McInnes } 11258965ea79SLois Curfman McInnes 11264a2ae208SSatish Balay #undef __FUNCT__ 11274a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 11285609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 11298965ea79SLois Curfman McInnes { 11308965ea79SLois Curfman McInnes Mat mat; 11313501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 113239ddd567SLois Curfman McInnes int ierr; 11338965ea79SLois Curfman McInnes 11343a40ed3dSBarry Smith PetscFunctionBegin; 11358965ea79SLois Curfman McInnes *newmat = 0; 1136273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr); 1137273d9f13SBarry Smith ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr); 1138b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1139b0a32e0cSBarry Smith mat->data = (void*)a; 1140549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 11413501a2bdSLois Curfman McInnes mat->factor = A->factor; 1142c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1143273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 11448965ea79SLois Curfman McInnes 11458965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 11468965ea79SLois Curfman McInnes a->rend = oldmat->rend; 11478965ea79SLois Curfman McInnes a->size = oldmat->size; 11488965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1149e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1150b9b97703SBarry Smith a->nvec = oldmat->nvec; 11513782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1152b0a32e0cSBarry Smith ierr = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1153b0a32e0cSBarry Smith PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1154549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 11558798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 11568965ea79SLois Curfman McInnes 1157329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 11585609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 1159b0a32e0cSBarry Smith PetscLogObjectParent(mat,a->A); 11608965ea79SLois Curfman McInnes *newmat = mat; 11613a40ed3dSBarry Smith PetscFunctionReturn(0); 11628965ea79SLois Curfman McInnes } 11638965ea79SLois Curfman McInnes 1164e090d566SSatish Balay #include "petscsys.h" 11658965ea79SLois Curfman McInnes 11664a2ae208SSatish Balay #undef __FUNCT__ 11674a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 116890ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat) 116990ace30eSBarry Smith { 117040011551SBarry Smith int *rowners,i,size,rank,m,ierr,nz,j; 1171*87828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 117290ace30eSBarry Smith MPI_Status status; 117390ace30eSBarry Smith 11743a40ed3dSBarry Smith PetscFunctionBegin; 1175d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1176d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 117790ace30eSBarry Smith 117890ace30eSBarry Smith /* determine ownership of all rows */ 117990ace30eSBarry Smith m = M/size + ((M % size) > rank); 1180b0a32e0cSBarry Smith ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1181ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 118290ace30eSBarry Smith rowners[0] = 0; 118390ace30eSBarry Smith for (i=2; i<=size; i++) { 118490ace30eSBarry Smith rowners[i] += rowners[i-1]; 118590ace30eSBarry Smith } 118690ace30eSBarry Smith 118790ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 118890ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 118990ace30eSBarry Smith 119090ace30eSBarry Smith if (!rank) { 1191*87828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 119290ace30eSBarry Smith 119390ace30eSBarry Smith /* read in my part of the matrix numerical values */ 11940752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 119590ace30eSBarry Smith 119690ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 119790ace30eSBarry Smith vals_ptr = vals; 119890ace30eSBarry Smith for (i=0; i<m; i++) { 119990ace30eSBarry Smith for (j=0; j<N; j++) { 120090ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 120190ace30eSBarry Smith } 120290ace30eSBarry Smith } 120390ace30eSBarry Smith 120490ace30eSBarry Smith /* read in other processors and ship out */ 120590ace30eSBarry Smith for (i=1; i<size; i++) { 120690ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 12070752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1208ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 120990ace30eSBarry Smith } 12103a40ed3dSBarry Smith } else { 121190ace30eSBarry Smith /* receive numeric values */ 1212*87828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 121390ace30eSBarry Smith 121490ace30eSBarry Smith /* receive message of values*/ 1215ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 121690ace30eSBarry Smith 121790ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 121890ace30eSBarry Smith vals_ptr = vals; 121990ace30eSBarry Smith for (i=0; i<m; i++) { 122090ace30eSBarry Smith for (j=0; j<N; j++) { 122190ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 122290ace30eSBarry Smith } 122390ace30eSBarry Smith } 122490ace30eSBarry Smith } 1225606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1226606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 12276d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12286d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12293a40ed3dSBarry Smith PetscFunctionReturn(0); 123090ace30eSBarry Smith } 123190ace30eSBarry Smith 1232273d9f13SBarry Smith EXTERN_C_BEGIN 12334a2ae208SSatish Balay #undef __FUNCT__ 12344a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1235b0a32e0cSBarry Smith int MatLoad_MPIDense(PetscViewer viewer,MatType type,Mat *newmat) 12368965ea79SLois Curfman McInnes { 12378965ea79SLois Curfman McInnes Mat A; 1238*87828ca2SBarry Smith PetscScalar *vals,*svals; 123919bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 12408965ea79SLois Curfman McInnes MPI_Status status; 12418965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 12428965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols; 124319bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 12443a40ed3dSBarry Smith int i,nz,ierr,j,rstart,rend,fd; 12458965ea79SLois Curfman McInnes 12463a40ed3dSBarry Smith PetscFunctionBegin; 1247d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1248d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 12498965ea79SLois Curfman McInnes if (!rank) { 1250b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 12510752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 125229bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 12538965ea79SLois Curfman McInnes } 12548965ea79SLois Curfman McInnes 1255ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 125690ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 125790ace30eSBarry Smith 125890ace30eSBarry Smith /* 125990ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 126090ace30eSBarry Smith */ 126190ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 12623a40ed3dSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 12633a40ed3dSBarry Smith PetscFunctionReturn(0); 126490ace30eSBarry Smith } 126590ace30eSBarry Smith 12668965ea79SLois Curfman McInnes /* determine ownership of all rows */ 12678965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 1268b0a32e0cSBarry Smith ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1269ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 12708965ea79SLois Curfman McInnes rowners[0] = 0; 12718965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 12728965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 12738965ea79SLois Curfman McInnes } 12748965ea79SLois Curfman McInnes rstart = rowners[rank]; 12758965ea79SLois Curfman McInnes rend = rowners[rank+1]; 12768965ea79SLois Curfman McInnes 12778965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 1278b0a32e0cSBarry Smith ierr = PetscMalloc(2*(rend-rstart)*sizeof(int),&ourlens);CHKERRQ(ierr); 12798965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 12808965ea79SLois Curfman McInnes if (!rank) { 1281b0a32e0cSBarry Smith ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr); 12820752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 1283b0a32e0cSBarry Smith ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr); 12848965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 1285ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1286606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1287ca161407SBarry Smith } else { 1288ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 12898965ea79SLois Curfman McInnes } 12908965ea79SLois Curfman McInnes 12918965ea79SLois Curfman McInnes if (!rank) { 12928965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 1293b0a32e0cSBarry Smith ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr); 1294549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 12958965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 12968965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 12978965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 12988965ea79SLois Curfman McInnes } 12998965ea79SLois Curfman McInnes } 1300606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 13018965ea79SLois Curfman McInnes 13028965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 13038965ea79SLois Curfman McInnes maxnz = 0; 13048965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 13050452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 13068965ea79SLois Curfman McInnes } 1307b0a32e0cSBarry Smith ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr); 13088965ea79SLois Curfman McInnes 13098965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 13108965ea79SLois Curfman McInnes nz = procsnz[0]; 1311b0a32e0cSBarry Smith ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr); 13120752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 13138965ea79SLois Curfman McInnes 13148965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 13158965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 13168965ea79SLois Curfman McInnes nz = procsnz[i]; 13170752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1318ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 13198965ea79SLois Curfman McInnes } 1320606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 13213a40ed3dSBarry Smith } else { 13228965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 13238965ea79SLois Curfman McInnes nz = 0; 13248965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13258965ea79SLois Curfman McInnes nz += ourlens[i]; 13268965ea79SLois Curfman McInnes } 1327b0a32e0cSBarry Smith ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr); 13288965ea79SLois Curfman McInnes 13298965ea79SLois Curfman McInnes /* receive message of column indices*/ 1330ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1331ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 133229bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 13338965ea79SLois Curfman McInnes } 13348965ea79SLois Curfman McInnes 13358965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1336549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 13378965ea79SLois Curfman McInnes jj = 0; 13388965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13398965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 13408965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 13418965ea79SLois Curfman McInnes jj++; 13428965ea79SLois Curfman McInnes } 13438965ea79SLois Curfman McInnes } 13448965ea79SLois Curfman McInnes 13458965ea79SLois Curfman McInnes /* create our matrix */ 13468965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13478965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 13488965ea79SLois Curfman McInnes } 1349b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 13508965ea79SLois Curfman McInnes A = *newmat; 13518965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13528965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 13538965ea79SLois Curfman McInnes } 13548965ea79SLois Curfman McInnes 13558965ea79SLois Curfman McInnes if (!rank) { 1356*87828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 13578965ea79SLois Curfman McInnes 13588965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 13598965ea79SLois Curfman McInnes nz = procsnz[0]; 13600752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 13618965ea79SLois Curfman McInnes 13628965ea79SLois Curfman McInnes /* insert into matrix */ 13638965ea79SLois Curfman McInnes jj = rstart; 13648965ea79SLois Curfman McInnes smycols = mycols; 13658965ea79SLois Curfman McInnes svals = vals; 13668965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13678965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 13688965ea79SLois Curfman McInnes smycols += ourlens[i]; 13698965ea79SLois Curfman McInnes svals += ourlens[i]; 13708965ea79SLois Curfman McInnes jj++; 13718965ea79SLois Curfman McInnes } 13728965ea79SLois Curfman McInnes 13738965ea79SLois Curfman McInnes /* read in other processors and ship out */ 13748965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 13758965ea79SLois Curfman McInnes nz = procsnz[i]; 13760752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1377ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 13788965ea79SLois Curfman McInnes } 1379606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 13803a40ed3dSBarry Smith } else { 13818965ea79SLois Curfman McInnes /* receive numeric values */ 1382*87828ca2SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 13838965ea79SLois Curfman McInnes 13848965ea79SLois Curfman McInnes /* receive message of values*/ 1385ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1386ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 138729bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 13888965ea79SLois Curfman McInnes 13898965ea79SLois Curfman McInnes /* insert into matrix */ 13908965ea79SLois Curfman McInnes jj = rstart; 13918965ea79SLois Curfman McInnes smycols = mycols; 13928965ea79SLois Curfman McInnes svals = vals; 13938965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13948965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 13958965ea79SLois Curfman McInnes smycols += ourlens[i]; 13968965ea79SLois Curfman McInnes svals += ourlens[i]; 13978965ea79SLois Curfman McInnes jj++; 13988965ea79SLois Curfman McInnes } 13998965ea79SLois Curfman McInnes } 1400606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1401606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1402606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1403606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 14048965ea79SLois Curfman McInnes 14056d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14066d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14073a40ed3dSBarry Smith PetscFunctionReturn(0); 14088965ea79SLois Curfman McInnes } 1409273d9f13SBarry Smith EXTERN_C_END 141090ace30eSBarry Smith 141190ace30eSBarry Smith 141290ace30eSBarry Smith 141390ace30eSBarry Smith 141490ace30eSBarry Smith 1415