1ed3cc1f0SBarry Smith /* 2ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 3ed3cc1f0SBarry Smith */ 4ed3cc1f0SBarry Smith 570f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 68965ea79SLois Curfman McInnes 70de54da6SSatish Balay EXTERN_C_BEGIN 84a2ae208SSatish Balay #undef __FUNCT__ 94a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 100de54da6SSatish Balay int MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 110de54da6SSatish Balay { 120de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 13cfce73b9SSatish Balay int m = A->m,rstart = mdn->rstart,ierr; 1487828ca2SBarry Smith PetscScalar *array; 150de54da6SSatish Balay MPI_Comm comm; 160de54da6SSatish Balay 170de54da6SSatish Balay PetscFunctionBegin; 18273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 190de54da6SSatish Balay 200de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 210de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 220de54da6SSatish Balay 230de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 240de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 255c5985e7SKris Buschelman ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr); 265c5985e7SKris Buschelman ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr); 275c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr); 280de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 290de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 300de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 310de54da6SSatish Balay 320de54da6SSatish Balay *iscopy = PETSC_TRUE; 330de54da6SSatish Balay PetscFunctionReturn(0); 340de54da6SSatish Balay } 350de54da6SSatish Balay EXTERN_C_END 360de54da6SSatish Balay 374a2ae208SSatish Balay #undef __FUNCT__ 384a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 39f15d580aSBarry Smith int MatSetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],const PetscScalar v[],InsertMode addv) 408965ea79SLois Curfman McInnes { 4139b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 4239b7565bSBarry Smith int ierr,i,j,rstart = A->rstart,rend = A->rend,row; 43273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 448965ea79SLois Curfman McInnes 453a40ed3dSBarry Smith PetscFunctionBegin; 468965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 475ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 48273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 498965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 508965ea79SLois Curfman McInnes row = idxm[i] - rstart; 5139b7565bSBarry Smith if (roworiented) { 5239b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 533a40ed3dSBarry Smith } else { 548965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 555ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 56273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 5739b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 5839b7565bSBarry Smith } 598965ea79SLois Curfman McInnes } 603a40ed3dSBarry Smith } else { 613782ba37SSatish Balay if (!A->donotstash) { 6239b7565bSBarry Smith if (roworiented) { 638798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 64d36fbae8SSatish Balay } else { 658798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 6639b7565bSBarry Smith } 67b49de8d1SLois Curfman McInnes } 68b49de8d1SLois Curfman McInnes } 693782ba37SSatish Balay } 703a40ed3dSBarry Smith PetscFunctionReturn(0); 71b49de8d1SLois Curfman McInnes } 72b49de8d1SLois Curfman McInnes 734a2ae208SSatish Balay #undef __FUNCT__ 744a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 75f15d580aSBarry Smith int MatGetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],PetscScalar v[]) 76b49de8d1SLois Curfman McInnes { 77b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 78b49de8d1SLois Curfman McInnes int ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row; 79b49de8d1SLois Curfman McInnes 803a40ed3dSBarry Smith PetscFunctionBegin; 81b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 8229bbc08cSBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 83273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 84b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 85b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 86b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 8729bbc08cSBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 88273d9f13SBarry Smith if (idxn[j] >= mat->N) { 8929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 90a8c6a408SBarry Smith } 91b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 92b49de8d1SLois Curfman McInnes } 93a8c6a408SBarry Smith } else { 9429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 958965ea79SLois Curfman McInnes } 968965ea79SLois Curfman McInnes } 973a40ed3dSBarry Smith PetscFunctionReturn(0); 988965ea79SLois Curfman McInnes } 998965ea79SLois Curfman McInnes 1004a2ae208SSatish Balay #undef __FUNCT__ 1014a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 1024e7234bfSBarry Smith int MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 103ff14e315SSatish Balay { 104ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 105ff14e315SSatish Balay int ierr; 106ff14e315SSatish Balay 1073a40ed3dSBarry Smith PetscFunctionBegin; 108ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1093a40ed3dSBarry Smith PetscFunctionReturn(0); 110ff14e315SSatish Balay } 111ff14e315SSatish Balay 1124a2ae208SSatish Balay #undef __FUNCT__ 1134a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 114ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 115ca3fa75bSLois Curfman McInnes { 116ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 117ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 118cfce73b9SSatish Balay int i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 11987828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 120ca3fa75bSLois Curfman McInnes Mat newmat; 121ca3fa75bSLois Curfman McInnes 122ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 123ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 124ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 125b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 126b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 127ca3fa75bSLois Curfman McInnes 128ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1297eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1307eba5e9cSLois Curfman McInnes original matrix! */ 131ca3fa75bSLois Curfman McInnes 132ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1337eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 134ca3fa75bSLois Curfman McInnes 135ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 136ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 13729bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 1387eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 139ca3fa75bSLois Curfman McInnes newmat = *B; 140ca3fa75bSLois Curfman McInnes } else { 141ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 142878740d9SKris Buschelman ierr = MatCreate(A->comm,nrows,cs,PETSC_DECIDE,ncols,&newmat);CHKERRQ(ierr); 143878740d9SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 144878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 145ca3fa75bSLois Curfman McInnes } 146ca3fa75bSLois Curfman McInnes 147ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 148ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 149ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 150ca3fa75bSLois Curfman McInnes 151ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 152ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 153ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 1547eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 155ca3fa75bSLois Curfman McInnes } 156ca3fa75bSLois Curfman McInnes } 157ca3fa75bSLois Curfman McInnes 158ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 159ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 160ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 161ca3fa75bSLois Curfman McInnes 162ca3fa75bSLois Curfman McInnes /* Free work space */ 163ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 164ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 165ca3fa75bSLois Curfman McInnes *B = newmat; 166ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 167ca3fa75bSLois Curfman McInnes } 168ca3fa75bSLois Curfman McInnes 1694a2ae208SSatish Balay #undef __FUNCT__ 1704a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 1714e7234bfSBarry Smith int MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 172ff14e315SSatish Balay { 1733a40ed3dSBarry Smith PetscFunctionBegin; 1743a40ed3dSBarry Smith PetscFunctionReturn(0); 175ff14e315SSatish Balay } 176ff14e315SSatish Balay 1774a2ae208SSatish Balay #undef __FUNCT__ 1784a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 1798f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 1808965ea79SLois Curfman McInnes { 18139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 1828965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 183d36fbae8SSatish Balay int ierr,nstash,reallocs; 1848965ea79SLois Curfman McInnes InsertMode addv; 1858965ea79SLois Curfman McInnes 1863a40ed3dSBarry Smith PetscFunctionBegin; 1878965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 188ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1897056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 19029bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 1918965ea79SLois Curfman McInnes } 192e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1938965ea79SLois Curfman McInnes 1948798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 1958798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 196b0a32e0cSBarry Smith PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs); 1973a40ed3dSBarry Smith PetscFunctionReturn(0); 1988965ea79SLois Curfman McInnes } 1998965ea79SLois Curfman McInnes 2004a2ae208SSatish Balay #undef __FUNCT__ 2014a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 2028f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2038965ea79SLois Curfman McInnes { 20439ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2057ef1d9bdSSatish Balay int i,n,ierr,*row,*col,flg,j,rstart,ncols; 20687828ca2SBarry Smith PetscScalar *val; 207e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2088965ea79SLois Curfman McInnes 2093a40ed3dSBarry Smith PetscFunctionBegin; 2108965ea79SLois Curfman McInnes /* wait on receives */ 2117ef1d9bdSSatish Balay while (1) { 2128798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2137ef1d9bdSSatish Balay if (!flg) break; 2148965ea79SLois Curfman McInnes 2157ef1d9bdSSatish Balay for (i=0; i<n;) { 2167ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2177ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2187ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2197ef1d9bdSSatish Balay else ncols = n-i; 2207ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2217ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2227ef1d9bdSSatish Balay i = j; 2238965ea79SLois Curfman McInnes } 2247ef1d9bdSSatish Balay } 2258798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2268965ea79SLois Curfman McInnes 22739ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 22839ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2298965ea79SLois Curfman McInnes 2306d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 23139ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2328965ea79SLois Curfman McInnes } 2333a40ed3dSBarry Smith PetscFunctionReturn(0); 2348965ea79SLois Curfman McInnes } 2358965ea79SLois Curfman McInnes 2364a2ae208SSatish Balay #undef __FUNCT__ 2374a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 2388f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A) 2398965ea79SLois Curfman McInnes { 2403a40ed3dSBarry Smith int ierr; 24139ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2423a40ed3dSBarry Smith 2433a40ed3dSBarry Smith PetscFunctionBegin; 2443a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2453a40ed3dSBarry Smith PetscFunctionReturn(0); 2468965ea79SLois Curfman McInnes } 2478965ea79SLois Curfman McInnes 2484a2ae208SSatish Balay #undef __FUNCT__ 2494a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIDense" 2508f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs) 2514e220ebcSLois Curfman McInnes { 2523a40ed3dSBarry Smith PetscFunctionBegin; 2534e220ebcSLois Curfman McInnes *bs = 1; 2543a40ed3dSBarry Smith PetscFunctionReturn(0); 2554e220ebcSLois Curfman McInnes } 2564e220ebcSLois Curfman McInnes 2578965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2588965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2598965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2603501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2618965ea79SLois Curfman McInnes routine. 2628965ea79SLois Curfman McInnes */ 2634a2ae208SSatish Balay #undef __FUNCT__ 2644a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 265268466fbSBarry Smith int MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag) 2668965ea79SLois Curfman McInnes { 26739ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2688965ea79SLois Curfman McInnes int i,ierr,N,*rows,*owners = l->rowners,size = l->size; 269c1dc657dSBarry Smith int *nprocs,j,idx,nsends; 2708965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2718965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2728965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2738965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2748965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2758965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2768965ea79SLois Curfman McInnes IS istmp; 27735d8aa7fSBarry Smith PetscTruth found; 2788965ea79SLois Curfman McInnes 2793a40ed3dSBarry Smith PetscFunctionBegin; 280b9b97703SBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 2818965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 2828965ea79SLois Curfman McInnes 2838965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 284b0a32e0cSBarry Smith ierr = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr); 285549d3d68SSatish Balay ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 286b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/ 2878965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 2888965ea79SLois Curfman McInnes idx = rows[i]; 28935d8aa7fSBarry Smith found = PETSC_FALSE; 2908965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 2918965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 292c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 2938965ea79SLois Curfman McInnes } 2948965ea79SLois Curfman McInnes } 29529bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 2968965ea79SLois Curfman McInnes } 297c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 2988965ea79SLois Curfman McInnes 2998965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 300c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3018965ea79SLois Curfman McInnes 3028965ea79SLois Curfman McInnes /* post receives: */ 303b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr); 304b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3058965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 306ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3078965ea79SLois Curfman McInnes } 3088965ea79SLois Curfman McInnes 3098965ea79SLois Curfman McInnes /* do sends: 3108965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3118965ea79SLois Curfman McInnes the ith processor 3128965ea79SLois Curfman McInnes */ 313b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr); 314b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 315b0a32e0cSBarry Smith ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr); 3168965ea79SLois Curfman McInnes starts[0] = 0; 317c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3188965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3198965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3208965ea79SLois Curfman McInnes } 3218965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3228965ea79SLois Curfman McInnes 3238965ea79SLois Curfman McInnes starts[0] = 0; 324c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3258965ea79SLois Curfman McInnes count = 0; 3268965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 327c1dc657dSBarry Smith if (nprocs[2*i+1]) { 328c1dc657dSBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3298965ea79SLois Curfman McInnes } 3308965ea79SLois Curfman McInnes } 331606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3328965ea79SLois Curfman McInnes 3338965ea79SLois Curfman McInnes base = owners[rank]; 3348965ea79SLois Curfman McInnes 3358965ea79SLois Curfman McInnes /* wait on receives */ 336b0a32e0cSBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr); 3378965ea79SLois Curfman McInnes source = lens + nrecvs; 3388965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3398965ea79SLois Curfman McInnes while (count) { 340ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3418965ea79SLois Curfman McInnes /* unpack receives into our local space */ 342ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 3438965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3448965ea79SLois Curfman McInnes lens[imdex] = n; 3458965ea79SLois Curfman McInnes slen += n; 3468965ea79SLois Curfman McInnes count--; 3478965ea79SLois Curfman McInnes } 348606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3498965ea79SLois Curfman McInnes 3508965ea79SLois Curfman McInnes /* move the data into the send scatter */ 351b0a32e0cSBarry Smith ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr); 3528965ea79SLois Curfman McInnes count = 0; 3538965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 3548965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3558965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 3568965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3578965ea79SLois Curfman McInnes } 3588965ea79SLois Curfman McInnes } 359606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 360606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 361606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 362606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3638965ea79SLois Curfman McInnes 3648965ea79SLois Curfman McInnes /* actually zap the local rows */ 365029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 366b0a32e0cSBarry Smith PetscLogObjectParent(A,istmp); 367606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3688965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 3698965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 3708965ea79SLois Curfman McInnes 3718965ea79SLois Curfman McInnes /* wait on sends */ 3728965ea79SLois Curfman McInnes if (nsends) { 373b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 374ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 375606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 3768965ea79SLois Curfman McInnes } 377606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 378606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 3798965ea79SLois Curfman McInnes 3803a40ed3dSBarry Smith PetscFunctionReturn(0); 3818965ea79SLois Curfman McInnes } 3828965ea79SLois Curfman McInnes 3834a2ae208SSatish Balay #undef __FUNCT__ 3844a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 3858f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3868965ea79SLois Curfman McInnes { 38739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 3888965ea79SLois Curfman McInnes int ierr; 389c456f294SBarry Smith 3903a40ed3dSBarry Smith PetscFunctionBegin; 39143a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39243a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39344cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 3943a40ed3dSBarry Smith PetscFunctionReturn(0); 3958965ea79SLois Curfman McInnes } 3968965ea79SLois Curfman McInnes 3974a2ae208SSatish Balay #undef __FUNCT__ 3984a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 3998f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4008965ea79SLois Curfman McInnes { 40139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4028965ea79SLois Curfman McInnes int ierr; 403c456f294SBarry Smith 4043a40ed3dSBarry Smith PetscFunctionBegin; 40543a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 40643a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 40744cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4083a40ed3dSBarry Smith PetscFunctionReturn(0); 4098965ea79SLois Curfman McInnes } 4108965ea79SLois Curfman McInnes 4114a2ae208SSatish Balay #undef __FUNCT__ 4124a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 4137c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 414096963f5SLois Curfman McInnes { 415096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 416096963f5SLois Curfman McInnes int ierr; 41787828ca2SBarry Smith PetscScalar zero = 0.0; 418096963f5SLois Curfman McInnes 4193a40ed3dSBarry Smith PetscFunctionBegin; 4203501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 4217c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 422537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 423537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4243a40ed3dSBarry Smith PetscFunctionReturn(0); 425096963f5SLois Curfman McInnes } 426096963f5SLois Curfman McInnes 4274a2ae208SSatish Balay #undef __FUNCT__ 4284a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 4297c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 430096963f5SLois Curfman McInnes { 431096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 432096963f5SLois Curfman McInnes int ierr; 433096963f5SLois Curfman McInnes 4343a40ed3dSBarry Smith PetscFunctionBegin; 4353501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4367c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 437537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 438537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4393a40ed3dSBarry Smith PetscFunctionReturn(0); 440096963f5SLois Curfman McInnes } 441096963f5SLois Curfman McInnes 4424a2ae208SSatish Balay #undef __FUNCT__ 4434a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 4448f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 4458965ea79SLois Curfman McInnes { 44639ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 447096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 448273d9f13SBarry Smith int ierr,len,i,n,m = A->m,radd; 44987828ca2SBarry Smith PetscScalar *x,zero = 0.0; 450ed3cc1f0SBarry Smith 4513a40ed3dSBarry Smith PetscFunctionBegin; 452273d9f13SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 4531ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 454096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 455273d9f13SBarry Smith if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 456273d9f13SBarry Smith len = PetscMin(a->A->m,a->A->n); 4577ddc982cSLois Curfman McInnes radd = a->rstart*m; 45844cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 459096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 460096963f5SLois Curfman McInnes } 4611ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 4623a40ed3dSBarry Smith PetscFunctionReturn(0); 4638965ea79SLois Curfman McInnes } 4648965ea79SLois Curfman McInnes 4654a2ae208SSatish Balay #undef __FUNCT__ 4664a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 467e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 4688965ea79SLois Curfman McInnes { 4693501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4708965ea79SLois Curfman McInnes int ierr; 471ed3cc1f0SBarry Smith 4723a40ed3dSBarry Smith PetscFunctionBegin; 47394d884c6SBarry Smith 474aa482453SBarry Smith #if defined(PETSC_USE_LOG) 475b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N); 4768965ea79SLois Curfman McInnes #endif 4778798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 478606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 4793501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 4803501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4813501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 482622d7880SLois Curfman McInnes if (mdn->factor) { 483606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 484606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 485606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 486606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 487622d7880SLois Curfman McInnes } 488606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 4893a40ed3dSBarry Smith PetscFunctionReturn(0); 4908965ea79SLois Curfman McInnes } 49139ddd567SLois Curfman McInnes 4924a2ae208SSatish Balay #undef __FUNCT__ 4934a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 494b0a32e0cSBarry Smith static int MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 4958965ea79SLois Curfman McInnes { 49639ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4978965ea79SLois Curfman McInnes int ierr; 4987056b6fcSBarry Smith 4993a40ed3dSBarry Smith PetscFunctionBegin; 50039ddd567SLois Curfman McInnes if (mdn->size == 1) { 50139ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5028965ea79SLois Curfman McInnes } 50329bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5043a40ed3dSBarry Smith PetscFunctionReturn(0); 5058965ea79SLois Curfman McInnes } 5068965ea79SLois Curfman McInnes 5074a2ae208SSatish Balay #undef __FUNCT__ 5084a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 509b0a32e0cSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5108965ea79SLois Curfman McInnes { 51139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 512fb9695e5SSatish Balay int ierr,size = mdn->size,rank = mdn->rank; 513b0a32e0cSBarry Smith PetscViewerType vtype; 514*32077d6dSBarry Smith PetscTruth iascii,isdraw; 515b0a32e0cSBarry Smith PetscViewer sviewer; 516f3ef73ceSBarry Smith PetscViewerFormat format; 5178965ea79SLois Curfman McInnes 5183a40ed3dSBarry Smith PetscFunctionBegin; 519*32077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 520fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 521*32077d6dSBarry Smith if (iascii) { 522b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 523b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 524456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 5254e220ebcSLois Curfman McInnes MatInfo info; 526888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 527b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m, 5286831982aSBarry Smith (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr); 529b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5303501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5313a40ed3dSBarry Smith PetscFunctionReturn(0); 532fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 5333a40ed3dSBarry Smith PetscFunctionReturn(0); 5348965ea79SLois Curfman McInnes } 535f1af5d2fSBarry Smith } else if (isdraw) { 536b0a32e0cSBarry Smith PetscDraw draw; 537f1af5d2fSBarry Smith PetscTruth isnull; 538f1af5d2fSBarry Smith 539b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 540b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 541f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 542f1af5d2fSBarry Smith } 54377ed5343SBarry Smith 5448965ea79SLois Curfman McInnes if (size == 1) { 54539ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5463a40ed3dSBarry Smith } else { 5478965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5488965ea79SLois Curfman McInnes Mat A; 549273d9f13SBarry Smith int M = mat->M,N = mat->N,m,row,i,nz,*cols; 55087828ca2SBarry Smith PetscScalar *vals; 5518965ea79SLois Curfman McInnes 5528965ea79SLois Curfman McInnes if (!rank) { 553878740d9SKris Buschelman ierr = MatCreate(mat->comm,M,N,M,N,&A);CHKERRQ(ierr); 5543a40ed3dSBarry Smith } else { 555878740d9SKris Buschelman ierr = MatCreate(mat->comm,0,0,M,N,&A);CHKERRQ(ierr); 5568965ea79SLois Curfman McInnes } 557be5d1d56SKris Buschelman /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */ 558878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 559878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 560b0a32e0cSBarry Smith PetscLogObjectParent(mat,A); 5618965ea79SLois Curfman McInnes 56239ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 56339ddd567SLois Curfman McInnes but it's quick for now */ 564273d9f13SBarry Smith row = mdn->rstart; m = mdn->A->m; 5658965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 56639ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 56739ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 56839ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 56939ddd567SLois Curfman McInnes row++; 5708965ea79SLois Curfman McInnes } 5718965ea79SLois Curfman McInnes 5726d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5736d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 574b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 575b9b97703SBarry Smith if (!rank) { 5766831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 5778965ea79SLois Curfman McInnes } 578b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 579b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5808965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 5818965ea79SLois Curfman McInnes } 5823a40ed3dSBarry Smith PetscFunctionReturn(0); 5838965ea79SLois Curfman McInnes } 5848965ea79SLois Curfman McInnes 5854a2ae208SSatish Balay #undef __FUNCT__ 5864a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 587b0a32e0cSBarry Smith int MatView_MPIDense(Mat mat,PetscViewer viewer) 5888965ea79SLois Curfman McInnes { 58939ddd567SLois Curfman McInnes int ierr; 590*32077d6dSBarry Smith PetscTruth iascii,isbinary,isdraw,issocket; 5918965ea79SLois Curfman McInnes 592433994e6SBarry Smith PetscFunctionBegin; 5930f5bd95cSBarry Smith 594*32077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 595fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 596b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 597fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 5980f5bd95cSBarry Smith 599*32077d6dSBarry Smith if (iascii || issocket || isdraw) { 600f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 6010f5bd95cSBarry Smith } else if (isbinary) { 6023a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6035cd90555SBarry Smith } else { 60429bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6058965ea79SLois Curfman McInnes } 6063a40ed3dSBarry Smith PetscFunctionReturn(0); 6078965ea79SLois Curfman McInnes } 6088965ea79SLois Curfman McInnes 6094a2ae208SSatish Balay #undef __FUNCT__ 6104a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 6118f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6128965ea79SLois Curfman McInnes { 6133501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 6143501a2bdSLois Curfman McInnes Mat mdn = mat->A; 6154e220ebcSLois Curfman McInnes int ierr; 616329f5518SBarry Smith PetscReal isend[5],irecv[5]; 6178965ea79SLois Curfman McInnes 6183a40ed3dSBarry Smith PetscFunctionBegin; 619273d9f13SBarry Smith info->rows_global = (double)A->M; 620273d9f13SBarry Smith info->columns_global = (double)A->N; 621273d9f13SBarry Smith info->rows_local = (double)A->m; 622273d9f13SBarry Smith info->columns_local = (double)A->N; 6234e220ebcSLois Curfman McInnes info->block_size = 1.0; 6244e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6254e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6264e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6278965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6284e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6294e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6304e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6314e220ebcSLois Curfman McInnes info->memory = isend[3]; 6324e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6338965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 634d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 6354e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6364e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6374e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6384e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6394e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6408965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 641d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 6424e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6434e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6444e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6454e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6464e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6478965ea79SLois Curfman McInnes } 6484e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6494e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6504e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6513a40ed3dSBarry Smith PetscFunctionReturn(0); 6528965ea79SLois Curfman McInnes } 6538965ea79SLois Curfman McInnes 6544a2ae208SSatish Balay #undef __FUNCT__ 6554a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 6568f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 6578965ea79SLois Curfman McInnes { 65839ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 659273d9f13SBarry Smith int ierr; 6608965ea79SLois Curfman McInnes 6613a40ed3dSBarry Smith PetscFunctionBegin; 66212c028f9SKris Buschelman switch (op) { 66312c028f9SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 66412c028f9SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 66512c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 66612c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 66712c028f9SKris Buschelman case MAT_COLUMNS_SORTED: 66812c028f9SKris Buschelman case MAT_COLUMNS_UNSORTED: 669273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 67012c028f9SKris Buschelman break; 67112c028f9SKris Buschelman case MAT_ROW_ORIENTED: 672273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 673273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 67412c028f9SKris Buschelman break; 67512c028f9SKris Buschelman case MAT_ROWS_SORTED: 67612c028f9SKris Buschelman case MAT_ROWS_UNSORTED: 67712c028f9SKris Buschelman case MAT_YES_NEW_DIAGONALS: 67812c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 679b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 68012c028f9SKris Buschelman break; 68112c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 682273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 683273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 68412c028f9SKris Buschelman break; 68512c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 686273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 68712c028f9SKris Buschelman break; 68812c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 68929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 69077e54ba9SKris Buschelman case MAT_SYMMETRIC: 69177e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 6929a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 6939a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 6949a4540c5SBarry Smith case MAT_HERMITIAN: 6959a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 6969a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 6979a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 69877e54ba9SKris Buschelman break; 69912c028f9SKris Buschelman default: 70029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 7013a40ed3dSBarry Smith } 7023a40ed3dSBarry Smith PetscFunctionReturn(0); 7038965ea79SLois Curfman McInnes } 7048965ea79SLois Curfman McInnes 7054a2ae208SSatish Balay #undef __FUNCT__ 7064a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIDense" 70787828ca2SBarry Smith int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v) 7088965ea79SLois Curfman McInnes { 7093501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7103a40ed3dSBarry Smith int lrow,rstart = mat->rstart,rend = mat->rend,ierr; 7118965ea79SLois Curfman McInnes 7123a40ed3dSBarry Smith PetscFunctionBegin; 71329bbc08cSBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 7148965ea79SLois Curfman McInnes lrow = row - rstart; 7153a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 7163a40ed3dSBarry Smith PetscFunctionReturn(0); 7178965ea79SLois Curfman McInnes } 7188965ea79SLois Curfman McInnes 7194a2ae208SSatish Balay #undef __FUNCT__ 7204a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIDense" 72187828ca2SBarry Smith int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v) 7228965ea79SLois Curfman McInnes { 723606d414cSSatish Balay int ierr; 724606d414cSSatish Balay 7253a40ed3dSBarry Smith PetscFunctionBegin; 726606d414cSSatish Balay if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 727606d414cSSatish Balay if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 7283a40ed3dSBarry Smith PetscFunctionReturn(0); 7298965ea79SLois Curfman McInnes } 7308965ea79SLois Curfman McInnes 7314a2ae208SSatish Balay #undef __FUNCT__ 7324a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 7335b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7345b2fa520SLois Curfman McInnes { 7355b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7365b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 73787828ca2SBarry Smith PetscScalar *l,*r,x,*v; 738273d9f13SBarry Smith int ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 7395b2fa520SLois Curfman McInnes 7405b2fa520SLois Curfman McInnes PetscFunctionBegin; 74172d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7425b2fa520SLois Curfman McInnes if (ll) { 74372d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 74429bbc08cSBarry Smith if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 7451ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7465b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7475b2fa520SLois Curfman McInnes x = l[i]; 7485b2fa520SLois Curfman McInnes v = mat->v + i; 7495b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7505b2fa520SLois Curfman McInnes } 7511ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 752b0a32e0cSBarry Smith PetscLogFlops(n*m); 7535b2fa520SLois Curfman McInnes } 7545b2fa520SLois Curfman McInnes if (rr) { 75572d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 75629bbc08cSBarry Smith if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 7575b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7585b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7591ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 7605b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 7615b2fa520SLois Curfman McInnes x = r[i]; 7625b2fa520SLois Curfman McInnes v = mat->v + i*m; 7635b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 7645b2fa520SLois Curfman McInnes } 7651ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 766b0a32e0cSBarry Smith PetscLogFlops(n*m); 7675b2fa520SLois Curfman McInnes } 7685b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 7695b2fa520SLois Curfman McInnes } 7705b2fa520SLois Curfman McInnes 7714a2ae208SSatish Balay #undef __FUNCT__ 7724a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 773064f8208SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 774096963f5SLois Curfman McInnes { 7753501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7763501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 7773501a2bdSLois Curfman McInnes int ierr,i,j; 778329f5518SBarry Smith PetscReal sum = 0.0; 77987828ca2SBarry Smith PetscScalar *v = mat->v; 7803501a2bdSLois Curfman McInnes 7813a40ed3dSBarry Smith PetscFunctionBegin; 7823501a2bdSLois Curfman McInnes if (mdn->size == 1) { 783064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 7843501a2bdSLois Curfman McInnes } else { 7853501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 786273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 787aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 788329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 7893501a2bdSLois Curfman McInnes #else 7903501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 7913501a2bdSLois Curfman McInnes #endif 7923501a2bdSLois Curfman McInnes } 793064f8208SBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 794064f8208SBarry Smith *nrm = sqrt(*nrm); 795b0a32e0cSBarry Smith PetscLogFlops(2*mdn->A->n*mdn->A->m); 7963a40ed3dSBarry Smith } else if (type == NORM_1) { 797329f5518SBarry Smith PetscReal *tmp,*tmp2; 798b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 799273d9f13SBarry Smith tmp2 = tmp + A->N; 800273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 801064f8208SBarry Smith *nrm = 0.0; 8023501a2bdSLois Curfman McInnes v = mat->v; 803273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 804273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 80567e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8063501a2bdSLois Curfman McInnes } 8073501a2bdSLois Curfman McInnes } 808d7d1e502SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 809273d9f13SBarry Smith for (j=0; j<A->N; j++) { 810064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8113501a2bdSLois Curfman McInnes } 812606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 813b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 8143a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 815329f5518SBarry Smith PetscReal ntemp; 8163501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 817064f8208SBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 8183a40ed3dSBarry Smith } else { 81929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8203501a2bdSLois Curfman McInnes } 8213501a2bdSLois Curfman McInnes } 8223a40ed3dSBarry Smith PetscFunctionReturn(0); 8233501a2bdSLois Curfman McInnes } 8243501a2bdSLois Curfman McInnes 8254a2ae208SSatish Balay #undef __FUNCT__ 8264a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 8278f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 8283501a2bdSLois Curfman McInnes { 8293501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8303501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8313501a2bdSLois Curfman McInnes Mat B; 832273d9f13SBarry Smith int M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 8333501a2bdSLois Curfman McInnes int j,i,ierr; 83487828ca2SBarry Smith PetscScalar *v; 8353501a2bdSLois Curfman McInnes 8363a40ed3dSBarry Smith PetscFunctionBegin; 8377c922b88SBarry Smith if (!matout && M != N) { 83829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 8397056b6fcSBarry Smith } 840878740d9SKris Buschelman ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B);CHKERRQ(ierr); 841878740d9SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 842878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr); 8433501a2bdSLois Curfman McInnes 844273d9f13SBarry Smith m = a->A->m; n = a->A->n; v = Aloc->v; 845b0a32e0cSBarry Smith ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr); 8463501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 8473501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8483501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8493501a2bdSLois Curfman McInnes v += m; 8503501a2bdSLois Curfman McInnes } 851606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 8526d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8536d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8547c922b88SBarry Smith if (matout) { 8553501a2bdSLois Curfman McInnes *matout = B; 8563501a2bdSLois Curfman McInnes } else { 857273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 8583501a2bdSLois Curfman McInnes } 8593a40ed3dSBarry Smith PetscFunctionReturn(0); 860096963f5SLois Curfman McInnes } 861096963f5SLois Curfman McInnes 862d9eff348SSatish Balay #include "petscblaslapack.h" 8634a2ae208SSatish Balay #undef __FUNCT__ 8644a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 865268466fbSBarry Smith int MatScale_MPIDense(const PetscScalar *alpha,Mat inA) 86644cd7ae7SLois Curfman McInnes { 86744cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 86844cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 86944cd7ae7SLois Curfman McInnes int one = 1,nz; 87044cd7ae7SLois Curfman McInnes 8713a40ed3dSBarry Smith PetscFunctionBegin; 872273d9f13SBarry Smith nz = inA->m*inA->N; 873268466fbSBarry Smith BLscal_(&nz,(PetscScalar*)alpha,a->v,&one); 874b0a32e0cSBarry Smith PetscLogFlops(nz); 8753a40ed3dSBarry Smith PetscFunctionReturn(0); 87644cd7ae7SLois Curfman McInnes } 87744cd7ae7SLois Curfman McInnes 8785609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 8798965ea79SLois Curfman McInnes 8804a2ae208SSatish Balay #undef __FUNCT__ 8814a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 882273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A) 883273d9f13SBarry Smith { 884273d9f13SBarry Smith int ierr; 885273d9f13SBarry Smith 886273d9f13SBarry Smith PetscFunctionBegin; 887273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 888273d9f13SBarry Smith PetscFunctionReturn(0); 889273d9f13SBarry Smith } 890273d9f13SBarry Smith 8918965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 89209dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 89309dc0095SBarry Smith MatGetRow_MPIDense, 89409dc0095SBarry Smith MatRestoreRow_MPIDense, 89509dc0095SBarry Smith MatMult_MPIDense, 89697304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 8977c922b88SBarry Smith MatMultTranspose_MPIDense, 8987c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 8998965ea79SLois Curfman McInnes 0, 90009dc0095SBarry Smith 0, 90109dc0095SBarry Smith 0, 90297304618SKris Buschelman /*10*/ 0, 90309dc0095SBarry Smith 0, 90409dc0095SBarry Smith 0, 90509dc0095SBarry Smith 0, 90609dc0095SBarry Smith MatTranspose_MPIDense, 90797304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 90897304618SKris Buschelman 0, 90909dc0095SBarry Smith MatGetDiagonal_MPIDense, 9105b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 91109dc0095SBarry Smith MatNorm_MPIDense, 91297304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 91309dc0095SBarry Smith MatAssemblyEnd_MPIDense, 91409dc0095SBarry Smith 0, 91509dc0095SBarry Smith MatSetOption_MPIDense, 91609dc0095SBarry Smith MatZeroEntries_MPIDense, 91797304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 91809dc0095SBarry Smith 0, 91909dc0095SBarry Smith 0, 92009dc0095SBarry Smith 0, 92109dc0095SBarry Smith 0, 92297304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 923273d9f13SBarry Smith 0, 92409dc0095SBarry Smith 0, 92509dc0095SBarry Smith MatGetArray_MPIDense, 92609dc0095SBarry Smith MatRestoreArray_MPIDense, 92797304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 92809dc0095SBarry Smith 0, 92909dc0095SBarry Smith 0, 93009dc0095SBarry Smith 0, 93109dc0095SBarry Smith 0, 93297304618SKris Buschelman /*40*/ 0, 9332ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 93409dc0095SBarry Smith 0, 93509dc0095SBarry Smith MatGetValues_MPIDense, 93609dc0095SBarry Smith 0, 93797304618SKris Buschelman /*45*/ 0, 93809dc0095SBarry Smith MatScale_MPIDense, 93909dc0095SBarry Smith 0, 94009dc0095SBarry Smith 0, 94109dc0095SBarry Smith 0, 94297304618SKris Buschelman /*50*/ MatGetBlockSize_MPIDense, 94309dc0095SBarry Smith 0, 94409dc0095SBarry Smith 0, 94509dc0095SBarry Smith 0, 94609dc0095SBarry Smith 0, 94797304618SKris Buschelman /*55*/ 0, 94809dc0095SBarry Smith 0, 94909dc0095SBarry Smith 0, 95009dc0095SBarry Smith 0, 95109dc0095SBarry Smith 0, 95297304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 953b9b97703SBarry Smith MatDestroy_MPIDense, 954b9b97703SBarry Smith MatView_MPIDense, 95597304618SKris Buschelman MatGetPetscMaps_Petsc, 95697304618SKris Buschelman 0, 95797304618SKris Buschelman /*65*/ 0, 95897304618SKris Buschelman 0, 95997304618SKris Buschelman 0, 96097304618SKris Buschelman 0, 96197304618SKris Buschelman 0, 96297304618SKris Buschelman /*70*/ 0, 96397304618SKris Buschelman 0, 96497304618SKris Buschelman 0, 96597304618SKris Buschelman 0, 96697304618SKris Buschelman 0, 96797304618SKris Buschelman /*75*/ 0, 96897304618SKris Buschelman 0, 96997304618SKris Buschelman 0, 97097304618SKris Buschelman 0, 97197304618SKris Buschelman 0, 97297304618SKris Buschelman /*80*/ 0, 97397304618SKris Buschelman 0, 97497304618SKris Buschelman 0, 97597304618SKris Buschelman 0, 97697304618SKris Buschelman /*85*/ MatLoad_MPIDense}; 9778965ea79SLois Curfman McInnes 978273d9f13SBarry Smith EXTERN_C_BEGIN 9794a2ae208SSatish Balay #undef __FUNCT__ 980a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 981a23d5eceSKris Buschelman int MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 982a23d5eceSKris Buschelman { 983a23d5eceSKris Buschelman Mat_MPIDense *a; 984a23d5eceSKris Buschelman int ierr; 985a23d5eceSKris Buschelman 986a23d5eceSKris Buschelman PetscFunctionBegin; 987a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 988a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 989a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 990a23d5eceSKris Buschelman 991a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 9925c5985e7SKris Buschelman ierr = MatCreate(PETSC_COMM_SELF,mat->m,mat->N,mat->m,mat->N,&a->A);CHKERRQ(ierr); 9935c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 9945c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 995a23d5eceSKris Buschelman PetscLogObjectParent(mat,a->A); 996a23d5eceSKris Buschelman PetscFunctionReturn(0); 997a23d5eceSKris Buschelman } 998a23d5eceSKris Buschelman EXTERN_C_END 999a23d5eceSKris Buschelman 10000bad9183SKris Buschelman /*MC 1001fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 10020bad9183SKris Buschelman 10030bad9183SKris Buschelman Options Database Keys: 10040bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 10050bad9183SKris Buschelman 10060bad9183SKris Buschelman Level: beginner 10070bad9183SKris Buschelman 10080bad9183SKris Buschelman .seealso: MatCreateMPIDense 10090bad9183SKris Buschelman M*/ 10100bad9183SKris Buschelman 1011a23d5eceSKris Buschelman EXTERN_C_BEGIN 1012a23d5eceSKris Buschelman #undef __FUNCT__ 10134a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1014273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat) 1015273d9f13SBarry Smith { 1016273d9f13SBarry Smith Mat_MPIDense *a; 1017273d9f13SBarry Smith int ierr,i; 1018273d9f13SBarry Smith 1019273d9f13SBarry Smith PetscFunctionBegin; 1020b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1021b0a32e0cSBarry Smith mat->data = (void*)a; 1022273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1023273d9f13SBarry Smith mat->factor = 0; 1024273d9f13SBarry Smith mat->mapping = 0; 1025273d9f13SBarry Smith 1026273d9f13SBarry Smith a->factor = 0; 1027273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1028273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1029273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1030273d9f13SBarry Smith 1031273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1032273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1033273d9f13SBarry Smith a->nvec = mat->n; 1034273d9f13SBarry Smith 1035273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 1036273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 10378a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 10388a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1039273d9f13SBarry Smith 1040273d9f13SBarry Smith /* build local table of row and column ownerships */ 1041b0a32e0cSBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1042273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 1043b0a32e0cSBarry Smith PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1044273d9f13SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1045273d9f13SBarry Smith a->rowners[0] = 0; 1046273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1047273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1048273d9f13SBarry Smith } 1049273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1050273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 1051273d9f13SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1052273d9f13SBarry Smith a->cowners[0] = 0; 1053273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1054273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1055273d9f13SBarry Smith } 1056273d9f13SBarry Smith 1057273d9f13SBarry Smith /* build cache for off array entries formed */ 1058273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1059273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1060273d9f13SBarry Smith 1061273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1062273d9f13SBarry Smith a->lvec = 0; 1063273d9f13SBarry Smith a->Mvctx = 0; 1064273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1065273d9f13SBarry Smith 1066273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1067273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1068273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1069a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1070a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1071a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1072273d9f13SBarry Smith PetscFunctionReturn(0); 1073273d9f13SBarry Smith } 1074273d9f13SBarry Smith EXTERN_C_END 1075273d9f13SBarry Smith 1076209238afSKris Buschelman /*MC 1077002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1078209238afSKris Buschelman 1079209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1080209238afSKris Buschelman and MATMPIDENSE otherwise. 1081209238afSKris Buschelman 1082209238afSKris Buschelman Options Database Keys: 1083209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1084209238afSKris Buschelman 1085209238afSKris Buschelman Level: beginner 1086209238afSKris Buschelman 1087209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1088209238afSKris Buschelman M*/ 1089209238afSKris Buschelman 1090209238afSKris Buschelman EXTERN_C_BEGIN 1091209238afSKris Buschelman #undef __FUNCT__ 1092209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1093209238afSKris Buschelman int MatCreate_Dense(Mat A) { 1094209238afSKris Buschelman int ierr,size; 1095209238afSKris Buschelman 1096209238afSKris Buschelman PetscFunctionBegin; 1097209238afSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1098209238afSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1099209238afSKris Buschelman if (size == 1) { 1100209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1101209238afSKris Buschelman } else { 1102209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1103209238afSKris Buschelman } 1104209238afSKris Buschelman PetscFunctionReturn(0); 1105209238afSKris Buschelman } 1106209238afSKris Buschelman EXTERN_C_END 1107209238afSKris Buschelman 11084a2ae208SSatish Balay #undef __FUNCT__ 11094a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1110273d9f13SBarry Smith /*@C 1111273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1112273d9f13SBarry Smith 1113273d9f13SBarry Smith Not collective 1114273d9f13SBarry Smith 1115273d9f13SBarry Smith Input Parameters: 1116273d9f13SBarry Smith . A - the matrix 1117273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1118273d9f13SBarry Smith to control all matrix memory allocation. 1119273d9f13SBarry Smith 1120273d9f13SBarry Smith Notes: 1121273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1122273d9f13SBarry Smith storage by columns. 1123273d9f13SBarry Smith 1124273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1125273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1126273d9f13SBarry Smith set data=PETSC_NULL. 1127273d9f13SBarry Smith 1128273d9f13SBarry Smith Level: intermediate 1129273d9f13SBarry Smith 1130273d9f13SBarry Smith .keywords: matrix,dense, parallel 1131273d9f13SBarry Smith 1132273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1133273d9f13SBarry Smith @*/ 113487828ca2SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1135273d9f13SBarry Smith { 1136a23d5eceSKris Buschelman int ierr,(*f)(Mat,PetscScalar *); 1137273d9f13SBarry Smith 1138273d9f13SBarry Smith PetscFunctionBegin; 1139565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1140a23d5eceSKris Buschelman if (f) { 1141a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1142a23d5eceSKris Buschelman } 1143273d9f13SBarry Smith PetscFunctionReturn(0); 1144273d9f13SBarry Smith } 1145273d9f13SBarry Smith 11464a2ae208SSatish Balay #undef __FUNCT__ 11474a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 11488965ea79SLois Curfman McInnes /*@C 114939ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 11508965ea79SLois Curfman McInnes 1151db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1152db81eaa0SLois Curfman McInnes 11538965ea79SLois Curfman McInnes Input Parameters: 1154db81eaa0SLois Curfman McInnes + comm - MPI communicator 11558965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1156db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 11578965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1158db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 11597f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1160dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 11618965ea79SLois Curfman McInnes 11628965ea79SLois Curfman McInnes Output Parameter: 1163477f1c0bSLois Curfman McInnes . A - the matrix 11648965ea79SLois Curfman McInnes 1165b259b22eSLois Curfman McInnes Notes: 116639ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 116739ddd567SLois Curfman McInnes storage by columns. 11688965ea79SLois Curfman McInnes 116918f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 117018f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 11717f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 117218f449edSLois Curfman McInnes 11738965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 11748965ea79SLois Curfman McInnes (possibly both). 11758965ea79SLois Curfman McInnes 1176027ccd11SLois Curfman McInnes Level: intermediate 1177027ccd11SLois Curfman McInnes 117839ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 11798965ea79SLois Curfman McInnes 118039ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 11818965ea79SLois Curfman McInnes @*/ 118287828ca2SBarry Smith int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A) 11838965ea79SLois Curfman McInnes { 1184273d9f13SBarry Smith int ierr,size; 11858965ea79SLois Curfman McInnes 11863a40ed3dSBarry Smith PetscFunctionBegin; 1187273d9f13SBarry Smith ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr); 1188273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1189273d9f13SBarry Smith if (size > 1) { 1190273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1191273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1192273d9f13SBarry Smith } else { 1193273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1194273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 11958c469469SLois Curfman McInnes } 11963a40ed3dSBarry Smith PetscFunctionReturn(0); 11978965ea79SLois Curfman McInnes } 11988965ea79SLois Curfman McInnes 11994a2ae208SSatish Balay #undef __FUNCT__ 12004a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 12015609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 12028965ea79SLois Curfman McInnes { 12038965ea79SLois Curfman McInnes Mat mat; 12043501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 120539ddd567SLois Curfman McInnes int ierr; 12068965ea79SLois Curfman McInnes 12073a40ed3dSBarry Smith PetscFunctionBegin; 12088965ea79SLois Curfman McInnes *newmat = 0; 1209273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr); 1210be5d1d56SKris Buschelman ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr); 1211b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1212b0a32e0cSBarry Smith mat->data = (void*)a; 1213549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 12143501a2bdSLois Curfman McInnes mat->factor = A->factor; 1215c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1216273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 12178965ea79SLois Curfman McInnes 12188965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 12198965ea79SLois Curfman McInnes a->rend = oldmat->rend; 12208965ea79SLois Curfman McInnes a->size = oldmat->size; 12218965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1222e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1223b9b97703SBarry Smith a->nvec = oldmat->nvec; 12243782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1225b0a32e0cSBarry Smith ierr = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1226b0a32e0cSBarry Smith PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1227549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 12288798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 12298965ea79SLois Curfman McInnes 1230329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 12315609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 1232b0a32e0cSBarry Smith PetscLogObjectParent(mat,a->A); 12338965ea79SLois Curfman McInnes *newmat = mat; 12343a40ed3dSBarry Smith PetscFunctionReturn(0); 12358965ea79SLois Curfman McInnes } 12368965ea79SLois Curfman McInnes 1237e090d566SSatish Balay #include "petscsys.h" 12388965ea79SLois Curfman McInnes 12394a2ae208SSatish Balay #undef __FUNCT__ 12404a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 1241be5d1d56SKris Buschelman int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,const MatType type,Mat *newmat) 124290ace30eSBarry Smith { 124340011551SBarry Smith int *rowners,i,size,rank,m,ierr,nz,j; 124487828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 124590ace30eSBarry Smith MPI_Status status; 124690ace30eSBarry Smith 12473a40ed3dSBarry Smith PetscFunctionBegin; 1248d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1249d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 125090ace30eSBarry Smith 125190ace30eSBarry Smith /* determine ownership of all rows */ 125290ace30eSBarry Smith m = M/size + ((M % size) > rank); 1253b0a32e0cSBarry Smith ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1254ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 125590ace30eSBarry Smith rowners[0] = 0; 125690ace30eSBarry Smith for (i=2; i<=size; i++) { 125790ace30eSBarry Smith rowners[i] += rowners[i-1]; 125890ace30eSBarry Smith } 125990ace30eSBarry Smith 1260878740d9SKris Buschelman ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 1261be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1262878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 126390ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 126490ace30eSBarry Smith 126590ace30eSBarry Smith if (!rank) { 126687828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 126790ace30eSBarry Smith 126890ace30eSBarry Smith /* read in my part of the matrix numerical values */ 12690752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 127090ace30eSBarry Smith 127190ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 127290ace30eSBarry Smith vals_ptr = vals; 127390ace30eSBarry Smith for (i=0; i<m; i++) { 127490ace30eSBarry Smith for (j=0; j<N; j++) { 127590ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 127690ace30eSBarry Smith } 127790ace30eSBarry Smith } 127890ace30eSBarry Smith 127990ace30eSBarry Smith /* read in other processors and ship out */ 128090ace30eSBarry Smith for (i=1; i<size; i++) { 128190ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 12820752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1283ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 128490ace30eSBarry Smith } 12853a40ed3dSBarry Smith } else { 128690ace30eSBarry Smith /* receive numeric values */ 128787828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 128890ace30eSBarry Smith 128990ace30eSBarry Smith /* receive message of values*/ 1290ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 129190ace30eSBarry Smith 129290ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 129390ace30eSBarry Smith vals_ptr = vals; 129490ace30eSBarry Smith for (i=0; i<m; i++) { 129590ace30eSBarry Smith for (j=0; j<N; j++) { 129690ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 129790ace30eSBarry Smith } 129890ace30eSBarry Smith } 129990ace30eSBarry Smith } 1300606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1301606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 13026d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13036d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13043a40ed3dSBarry Smith PetscFunctionReturn(0); 130590ace30eSBarry Smith } 130690ace30eSBarry Smith 13074a2ae208SSatish Balay #undef __FUNCT__ 13084a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 13098e9aea5cSBarry Smith int MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat) 13108965ea79SLois Curfman McInnes { 13118965ea79SLois Curfman McInnes Mat A; 131287828ca2SBarry Smith PetscScalar *vals,*svals; 131319bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 13148965ea79SLois Curfman McInnes MPI_Status status; 13158965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 13168965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols; 131719bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 13183a40ed3dSBarry Smith int i,nz,ierr,j,rstart,rend,fd; 13198965ea79SLois Curfman McInnes 13203a40ed3dSBarry Smith PetscFunctionBegin; 1321d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1322d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13238965ea79SLois Curfman McInnes if (!rank) { 1324b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 13250752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1326552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 13278965ea79SLois Curfman McInnes } 13288965ea79SLois Curfman McInnes 1329ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 133090ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 133190ace30eSBarry Smith 133290ace30eSBarry Smith /* 133390ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 133490ace30eSBarry Smith */ 133590ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1336be5d1d56SKris Buschelman ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 13373a40ed3dSBarry Smith PetscFunctionReturn(0); 133890ace30eSBarry Smith } 133990ace30eSBarry Smith 13408965ea79SLois Curfman McInnes /* determine ownership of all rows */ 13418965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 1342b0a32e0cSBarry Smith ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1343ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 13448965ea79SLois Curfman McInnes rowners[0] = 0; 13458965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 13468965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 13478965ea79SLois Curfman McInnes } 13488965ea79SLois Curfman McInnes rstart = rowners[rank]; 13498965ea79SLois Curfman McInnes rend = rowners[rank+1]; 13508965ea79SLois Curfman McInnes 13518965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 135284d528b1SBarry Smith ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(int),&ourlens);CHKERRQ(ierr); 13538965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 13548965ea79SLois Curfman McInnes if (!rank) { 1355b0a32e0cSBarry Smith ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr); 13560752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 1357b0a32e0cSBarry Smith ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr); 13588965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 1359ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1360606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1361ca161407SBarry Smith } else { 1362ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 13638965ea79SLois Curfman McInnes } 13648965ea79SLois Curfman McInnes 13658965ea79SLois Curfman McInnes if (!rank) { 13668965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 1367b0a32e0cSBarry Smith ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr); 1368549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 13698965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 13708965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 13718965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 13728965ea79SLois Curfman McInnes } 13738965ea79SLois Curfman McInnes } 1374606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 13758965ea79SLois Curfman McInnes 13768965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 13778965ea79SLois Curfman McInnes maxnz = 0; 13788965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 13790452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 13808965ea79SLois Curfman McInnes } 1381b0a32e0cSBarry Smith ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr); 13828965ea79SLois Curfman McInnes 13838965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 13848965ea79SLois Curfman McInnes nz = procsnz[0]; 1385b0a32e0cSBarry Smith ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr); 13860752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 13878965ea79SLois Curfman McInnes 13888965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 13898965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 13908965ea79SLois Curfman McInnes nz = procsnz[i]; 13910752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1392ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 13938965ea79SLois Curfman McInnes } 1394606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 13953a40ed3dSBarry Smith } else { 13968965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 13978965ea79SLois Curfman McInnes nz = 0; 13988965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13998965ea79SLois Curfman McInnes nz += ourlens[i]; 14008965ea79SLois Curfman McInnes } 140184d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(int),&mycols);CHKERRQ(ierr); 14028965ea79SLois Curfman McInnes 14038965ea79SLois Curfman McInnes /* receive message of column indices*/ 1404ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1405ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 140629bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14078965ea79SLois Curfman McInnes } 14088965ea79SLois Curfman McInnes 14098965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1410549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 14118965ea79SLois Curfman McInnes jj = 0; 14128965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14138965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 14148965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 14158965ea79SLois Curfman McInnes jj++; 14168965ea79SLois Curfman McInnes } 14178965ea79SLois Curfman McInnes } 14188965ea79SLois Curfman McInnes 14198965ea79SLois Curfman McInnes /* create our matrix */ 14208965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14218965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 14228965ea79SLois Curfman McInnes } 1423878740d9SKris Buschelman ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 1424be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1425878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 14268965ea79SLois Curfman McInnes A = *newmat; 14278965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14288965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 14298965ea79SLois Curfman McInnes } 14308965ea79SLois Curfman McInnes 14318965ea79SLois Curfman McInnes if (!rank) { 143287828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14338965ea79SLois Curfman McInnes 14348965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 14358965ea79SLois Curfman McInnes nz = procsnz[0]; 14360752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 14378965ea79SLois Curfman McInnes 14388965ea79SLois Curfman McInnes /* insert into matrix */ 14398965ea79SLois Curfman McInnes jj = rstart; 14408965ea79SLois Curfman McInnes smycols = mycols; 14418965ea79SLois Curfman McInnes svals = vals; 14428965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14438965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14448965ea79SLois Curfman McInnes smycols += ourlens[i]; 14458965ea79SLois Curfman McInnes svals += ourlens[i]; 14468965ea79SLois Curfman McInnes jj++; 14478965ea79SLois Curfman McInnes } 14488965ea79SLois Curfman McInnes 14498965ea79SLois Curfman McInnes /* read in other processors and ship out */ 14508965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14518965ea79SLois Curfman McInnes nz = procsnz[i]; 14520752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1453ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 14548965ea79SLois Curfman McInnes } 1455606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 14563a40ed3dSBarry Smith } else { 14578965ea79SLois Curfman McInnes /* receive numeric values */ 145884d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14598965ea79SLois Curfman McInnes 14608965ea79SLois Curfman McInnes /* receive message of values*/ 1461ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1462ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 146329bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14648965ea79SLois Curfman McInnes 14658965ea79SLois Curfman McInnes /* insert into matrix */ 14668965ea79SLois Curfman McInnes jj = rstart; 14678965ea79SLois Curfman McInnes smycols = mycols; 14688965ea79SLois Curfman McInnes svals = vals; 14698965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14708965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14718965ea79SLois Curfman McInnes smycols += ourlens[i]; 14728965ea79SLois Curfman McInnes svals += ourlens[i]; 14738965ea79SLois Curfman McInnes jj++; 14748965ea79SLois Curfman McInnes } 14758965ea79SLois Curfman McInnes } 1476606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1477606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1478606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1479606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 14808965ea79SLois Curfman McInnes 14816d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14826d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14833a40ed3dSBarry Smith PetscFunctionReturn(0); 14848965ea79SLois Curfman McInnes } 148590ace30eSBarry Smith 148690ace30eSBarry Smith 148790ace30eSBarry Smith 148890ace30eSBarry Smith 148990ace30eSBarry Smith 1490