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" 10dfbe8321SBarry Smith PetscErrorCode MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 110de54da6SSatish Balay { 120de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 136849ba73SBarry Smith PetscErrorCode ierr; 146849ba73SBarry Smith int m = A->m,rstart = mdn->rstart; 1587828ca2SBarry Smith PetscScalar *array; 160de54da6SSatish Balay MPI_Comm comm; 170de54da6SSatish Balay 180de54da6SSatish Balay PetscFunctionBegin; 19273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 200de54da6SSatish Balay 210de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 220de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 230de54da6SSatish Balay 240de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 250de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 265c5985e7SKris Buschelman ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr); 275c5985e7SKris Buschelman ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr); 285c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);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 384a2ae208SSatish Balay #undef __FUNCT__ 394a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 40dfbe8321SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],const PetscScalar v[],InsertMode addv) 418965ea79SLois Curfman McInnes { 4239b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 43dfbe8321SBarry Smith PetscErrorCode ierr; 44dfbe8321SBarry Smith int i,j,rstart = A->rstart,rend = A->rend,row; 45273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 468965ea79SLois Curfman McInnes 473a40ed3dSBarry Smith PetscFunctionBegin; 488965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 495ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 50273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 518965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 528965ea79SLois Curfman McInnes row = idxm[i] - rstart; 5339b7565bSBarry Smith if (roworiented) { 5439b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 553a40ed3dSBarry Smith } else { 568965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 575ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 58273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 5939b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 6039b7565bSBarry Smith } 618965ea79SLois Curfman McInnes } 623a40ed3dSBarry Smith } else { 633782ba37SSatish Balay if (!A->donotstash) { 6439b7565bSBarry Smith if (roworiented) { 658798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 66d36fbae8SSatish Balay } else { 678798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 6839b7565bSBarry Smith } 69b49de8d1SLois Curfman McInnes } 70b49de8d1SLois Curfman McInnes } 713782ba37SSatish Balay } 723a40ed3dSBarry Smith PetscFunctionReturn(0); 73b49de8d1SLois Curfman McInnes } 74b49de8d1SLois Curfman McInnes 754a2ae208SSatish Balay #undef __FUNCT__ 764a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 77dfbe8321SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],PetscScalar v[]) 78b49de8d1SLois Curfman McInnes { 79b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 80dfbe8321SBarry Smith PetscErrorCode ierr; 81dfbe8321SBarry Smith int 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" 105dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 106ff14e315SSatish Balay { 107ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 108dfbe8321SBarry Smith PetscErrorCode 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" 1176849ba73SBarry Smith static PetscErrorCode 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; 1216849ba73SBarry Smith PetscErrorCode ierr; 1226849ba73SBarry Smith int i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 12387828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 124ca3fa75bSLois Curfman McInnes Mat newmat; 125ca3fa75bSLois Curfman McInnes 126ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 127ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 128ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 129b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 130b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 131ca3fa75bSLois Curfman McInnes 132ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1337eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1347eba5e9cSLois Curfman McInnes original matrix! */ 135ca3fa75bSLois Curfman McInnes 136ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1377eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 138ca3fa75bSLois Curfman McInnes 139ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 140ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 14129bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 1427eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 143ca3fa75bSLois Curfman McInnes newmat = *B; 144ca3fa75bSLois Curfman McInnes } else { 145ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 146878740d9SKris Buschelman ierr = MatCreate(A->comm,nrows,cs,PETSC_DECIDE,ncols,&newmat);CHKERRQ(ierr); 147878740d9SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 148878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 149ca3fa75bSLois Curfman McInnes } 150ca3fa75bSLois Curfman McInnes 151ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 152ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 153ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 154ca3fa75bSLois Curfman McInnes 155ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 156ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 157ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 1587eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 159ca3fa75bSLois Curfman McInnes } 160ca3fa75bSLois Curfman McInnes } 161ca3fa75bSLois Curfman McInnes 162ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 163ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165ca3fa75bSLois Curfman McInnes 166ca3fa75bSLois Curfman McInnes /* Free work space */ 167ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 168ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 169ca3fa75bSLois Curfman McInnes *B = newmat; 170ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 171ca3fa75bSLois Curfman McInnes } 172ca3fa75bSLois Curfman McInnes 1734a2ae208SSatish Balay #undef __FUNCT__ 1744a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 175dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 176ff14e315SSatish Balay { 1773a40ed3dSBarry Smith PetscFunctionBegin; 1783a40ed3dSBarry Smith PetscFunctionReturn(0); 179ff14e315SSatish Balay } 180ff14e315SSatish Balay 1814a2ae208SSatish Balay #undef __FUNCT__ 1824a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 183dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 1848965ea79SLois Curfman McInnes { 18539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 1868965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 187dfbe8321SBarry Smith PetscErrorCode ierr; 188dfbe8321SBarry Smith int nstash,reallocs; 1898965ea79SLois Curfman McInnes InsertMode addv; 1908965ea79SLois Curfman McInnes 1913a40ed3dSBarry Smith PetscFunctionBegin; 1928965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 193ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1947056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 19529bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 1968965ea79SLois Curfman McInnes } 197e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1988965ea79SLois Curfman McInnes 1998798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 2008798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 20177431f27SBarry Smith PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs); 2023a40ed3dSBarry Smith PetscFunctionReturn(0); 2038965ea79SLois Curfman McInnes } 2048965ea79SLois Curfman McInnes 2054a2ae208SSatish Balay #undef __FUNCT__ 2064a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 207dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2088965ea79SLois Curfman McInnes { 20939ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2106849ba73SBarry Smith PetscErrorCode ierr; 2116849ba73SBarry Smith int i,n,*row,*col,flg,j,rstart,ncols; 21287828ca2SBarry Smith PetscScalar *val; 213e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2148965ea79SLois Curfman McInnes 2153a40ed3dSBarry Smith PetscFunctionBegin; 2168965ea79SLois Curfman McInnes /* wait on receives */ 2177ef1d9bdSSatish Balay while (1) { 2188798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2197ef1d9bdSSatish Balay if (!flg) break; 2208965ea79SLois Curfman McInnes 2217ef1d9bdSSatish Balay for (i=0; i<n;) { 2227ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2237ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2247ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2257ef1d9bdSSatish Balay else ncols = n-i; 2267ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2277ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2287ef1d9bdSSatish Balay i = j; 2298965ea79SLois Curfman McInnes } 2307ef1d9bdSSatish Balay } 2318798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2328965ea79SLois Curfman McInnes 23339ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 23439ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2358965ea79SLois Curfman McInnes 2366d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 23739ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2388965ea79SLois Curfman McInnes } 2393a40ed3dSBarry Smith PetscFunctionReturn(0); 2408965ea79SLois Curfman McInnes } 2418965ea79SLois Curfman McInnes 2424a2ae208SSatish Balay #undef __FUNCT__ 2434a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 244dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 2458965ea79SLois Curfman McInnes { 246dfbe8321SBarry Smith PetscErrorCode ierr; 24739ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2483a40ed3dSBarry Smith 2493a40ed3dSBarry Smith PetscFunctionBegin; 2503a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2513a40ed3dSBarry Smith PetscFunctionReturn(0); 2528965ea79SLois Curfman McInnes } 2538965ea79SLois Curfman McInnes 2548965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2558965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2568965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2573501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2588965ea79SLois Curfman McInnes routine. 2598965ea79SLois Curfman McInnes */ 2604a2ae208SSatish Balay #undef __FUNCT__ 2614a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 262dfbe8321SBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag) 2638965ea79SLois Curfman McInnes { 26439ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2656849ba73SBarry Smith PetscErrorCode ierr; 2666849ba73SBarry Smith int i,N,*rows,*owners = l->rowners,size = l->size; 267c1dc657dSBarry Smith int *nprocs,j,idx,nsends; 2688965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2698965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2708965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2718965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2728965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2738965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2748965ea79SLois Curfman McInnes IS istmp; 27535d8aa7fSBarry Smith PetscTruth found; 2768965ea79SLois Curfman McInnes 2773a40ed3dSBarry Smith PetscFunctionBegin; 278b9b97703SBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 2798965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 2808965ea79SLois Curfman McInnes 2818965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 282b0a32e0cSBarry Smith ierr = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr); 283549d3d68SSatish Balay ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 284b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/ 2858965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 2868965ea79SLois Curfman McInnes idx = rows[i]; 28735d8aa7fSBarry Smith found = PETSC_FALSE; 2888965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 2898965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 290c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 2918965ea79SLois Curfman McInnes } 2928965ea79SLois Curfman McInnes } 29329bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 2948965ea79SLois Curfman McInnes } 295c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 2968965ea79SLois Curfman McInnes 2978965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 298c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 2998965ea79SLois Curfman McInnes 3008965ea79SLois Curfman McInnes /* post receives: */ 301b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr); 302b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3038965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 304ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3058965ea79SLois Curfman McInnes } 3068965ea79SLois Curfman McInnes 3078965ea79SLois Curfman McInnes /* do sends: 3088965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3098965ea79SLois Curfman McInnes the ith processor 3108965ea79SLois Curfman McInnes */ 311b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr); 312b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 313b0a32e0cSBarry Smith ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr); 3148965ea79SLois Curfman McInnes starts[0] = 0; 315c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3168965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3178965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3188965ea79SLois Curfman McInnes } 3198965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3208965ea79SLois Curfman McInnes 3218965ea79SLois Curfman McInnes starts[0] = 0; 322c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3238965ea79SLois Curfman McInnes count = 0; 3248965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 325c1dc657dSBarry Smith if (nprocs[2*i+1]) { 326c1dc657dSBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3278965ea79SLois Curfman McInnes } 3288965ea79SLois Curfman McInnes } 329606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3308965ea79SLois Curfman McInnes 3318965ea79SLois Curfman McInnes base = owners[rank]; 3328965ea79SLois Curfman McInnes 3338965ea79SLois Curfman McInnes /* wait on receives */ 334b0a32e0cSBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr); 3358965ea79SLois Curfman McInnes source = lens + nrecvs; 3368965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3378965ea79SLois Curfman McInnes while (count) { 338ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3398965ea79SLois Curfman McInnes /* unpack receives into our local space */ 340ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 3418965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3428965ea79SLois Curfman McInnes lens[imdex] = n; 3438965ea79SLois Curfman McInnes slen += n; 3448965ea79SLois Curfman McInnes count--; 3458965ea79SLois Curfman McInnes } 346606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3478965ea79SLois Curfman McInnes 3488965ea79SLois Curfman McInnes /* move the data into the send scatter */ 349b0a32e0cSBarry Smith ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr); 3508965ea79SLois Curfman McInnes count = 0; 3518965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 3528965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3538965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 3548965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3558965ea79SLois Curfman McInnes } 3568965ea79SLois Curfman McInnes } 357606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 358606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 359606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 360606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3618965ea79SLois Curfman McInnes 3628965ea79SLois Curfman McInnes /* actually zap the local rows */ 363029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 364b0a32e0cSBarry Smith PetscLogObjectParent(A,istmp); 365606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3668965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 3678965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 3688965ea79SLois Curfman McInnes 3698965ea79SLois Curfman McInnes /* wait on sends */ 3708965ea79SLois Curfman McInnes if (nsends) { 371b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 372ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 373606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 3748965ea79SLois Curfman McInnes } 375606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 376606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 3778965ea79SLois Curfman McInnes 3783a40ed3dSBarry Smith PetscFunctionReturn(0); 3798965ea79SLois Curfman McInnes } 3808965ea79SLois Curfman McInnes 3814a2ae208SSatish Balay #undef __FUNCT__ 3824a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 383dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3848965ea79SLois Curfman McInnes { 38539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 386dfbe8321SBarry Smith PetscErrorCode ierr; 387c456f294SBarry Smith 3883a40ed3dSBarry Smith PetscFunctionBegin; 38943a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39043a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39144cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 3923a40ed3dSBarry Smith PetscFunctionReturn(0); 3938965ea79SLois Curfman McInnes } 3948965ea79SLois Curfman McInnes 3954a2ae208SSatish Balay #undef __FUNCT__ 3964a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 397dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3988965ea79SLois Curfman McInnes { 39939ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 400dfbe8321SBarry Smith PetscErrorCode ierr; 401c456f294SBarry Smith 4023a40ed3dSBarry Smith PetscFunctionBegin; 40343a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 40443a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 40544cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4063a40ed3dSBarry Smith PetscFunctionReturn(0); 4078965ea79SLois Curfman McInnes } 4088965ea79SLois Curfman McInnes 4094a2ae208SSatish Balay #undef __FUNCT__ 4104a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 411dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 412096963f5SLois Curfman McInnes { 413096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 414dfbe8321SBarry Smith PetscErrorCode ierr; 41587828ca2SBarry Smith PetscScalar zero = 0.0; 416096963f5SLois Curfman McInnes 4173a40ed3dSBarry Smith PetscFunctionBegin; 4183501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 4197c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 420537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 421537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4223a40ed3dSBarry Smith PetscFunctionReturn(0); 423096963f5SLois Curfman McInnes } 424096963f5SLois Curfman McInnes 4254a2ae208SSatish Balay #undef __FUNCT__ 4264a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 427dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 428096963f5SLois Curfman McInnes { 429096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 430dfbe8321SBarry Smith PetscErrorCode ierr; 431096963f5SLois Curfman McInnes 4323a40ed3dSBarry Smith PetscFunctionBegin; 4333501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4347c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 435537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 436537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4373a40ed3dSBarry Smith PetscFunctionReturn(0); 438096963f5SLois Curfman McInnes } 439096963f5SLois Curfman McInnes 4404a2ae208SSatish Balay #undef __FUNCT__ 4414a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 442dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 4438965ea79SLois Curfman McInnes { 44439ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 445096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 446dfbe8321SBarry Smith PetscErrorCode ierr; 447dfbe8321SBarry Smith int len,i,n,m = A->m,radd; 44887828ca2SBarry Smith PetscScalar *x,zero = 0.0; 449ed3cc1f0SBarry Smith 4503a40ed3dSBarry Smith PetscFunctionBegin; 451273d9f13SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 4521ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 453096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 454273d9f13SBarry Smith if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 455273d9f13SBarry Smith len = PetscMin(a->A->m,a->A->n); 4567ddc982cSLois Curfman McInnes radd = a->rstart*m; 45744cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 458096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 459096963f5SLois Curfman McInnes } 4601ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 4613a40ed3dSBarry Smith PetscFunctionReturn(0); 4628965ea79SLois Curfman McInnes } 4638965ea79SLois Curfman McInnes 4644a2ae208SSatish Balay #undef __FUNCT__ 4654a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 466dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 4678965ea79SLois Curfman McInnes { 4683501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 469dfbe8321SBarry Smith PetscErrorCode ierr; 470ed3cc1f0SBarry Smith 4713a40ed3dSBarry Smith PetscFunctionBegin; 47294d884c6SBarry Smith 473aa482453SBarry Smith #if defined(PETSC_USE_LOG) 47477431f27SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N); 4758965ea79SLois Curfman McInnes #endif 4768798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 477606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 4783501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 4793501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4803501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 481622d7880SLois Curfman McInnes if (mdn->factor) { 482606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 483606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 484606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 485606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 486622d7880SLois Curfman McInnes } 487606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 488901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr); 489901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 4903a40ed3dSBarry Smith PetscFunctionReturn(0); 4918965ea79SLois Curfman McInnes } 49239ddd567SLois Curfman McInnes 4934a2ae208SSatish Balay #undef __FUNCT__ 4944a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 4956849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 4968965ea79SLois Curfman McInnes { 49739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 498dfbe8321SBarry Smith PetscErrorCode ierr; 4997056b6fcSBarry Smith 5003a40ed3dSBarry Smith PetscFunctionBegin; 50139ddd567SLois Curfman McInnes if (mdn->size == 1) { 50239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5038965ea79SLois Curfman McInnes } 50429bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5053a40ed3dSBarry Smith PetscFunctionReturn(0); 5068965ea79SLois Curfman McInnes } 5078965ea79SLois Curfman McInnes 5084a2ae208SSatish Balay #undef __FUNCT__ 5094a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 5106849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5118965ea79SLois Curfman McInnes { 51239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 513dfbe8321SBarry Smith PetscErrorCode ierr; 51432dcc486SBarry Smith PetscMPIInt size = mdn->size,rank = mdn->rank; 515b0a32e0cSBarry Smith PetscViewerType vtype; 51632077d6dSBarry Smith PetscTruth iascii,isdraw; 517b0a32e0cSBarry Smith PetscViewer sviewer; 518f3ef73ceSBarry Smith PetscViewerFormat format; 5198965ea79SLois Curfman McInnes 5203a40ed3dSBarry Smith PetscFunctionBegin; 52132077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 522fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 52332077d6dSBarry Smith if (iascii) { 524b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 525b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 526456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 5274e220ebcSLois Curfman McInnes MatInfo info; 528888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 52977431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m, 53077431f27SBarry Smith (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 531b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5323501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5333a40ed3dSBarry Smith PetscFunctionReturn(0); 534fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 5353a40ed3dSBarry Smith PetscFunctionReturn(0); 5368965ea79SLois Curfman McInnes } 537f1af5d2fSBarry Smith } else if (isdraw) { 538b0a32e0cSBarry Smith PetscDraw draw; 539f1af5d2fSBarry Smith PetscTruth isnull; 540f1af5d2fSBarry Smith 541b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 542b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 543f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 544f1af5d2fSBarry Smith } 54577ed5343SBarry Smith 5468965ea79SLois Curfman McInnes if (size == 1) { 54739ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5483a40ed3dSBarry Smith } else { 5498965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5508965ea79SLois Curfman McInnes Mat A; 551b3cc6726SBarry Smith int M = mat->M,N = mat->N,m,row,i,nz; 552b3cc6726SBarry Smith const int *cols; 553b3cc6726SBarry Smith const PetscScalar *vals; 5548965ea79SLois Curfman McInnes 5558965ea79SLois Curfman McInnes if (!rank) { 556878740d9SKris Buschelman ierr = MatCreate(mat->comm,M,N,M,N,&A);CHKERRQ(ierr); 5573a40ed3dSBarry Smith } else { 558878740d9SKris Buschelman ierr = MatCreate(mat->comm,0,0,M,N,&A);CHKERRQ(ierr); 5598965ea79SLois Curfman McInnes } 560be5d1d56SKris Buschelman /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */ 561878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 562878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 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" 590dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 5918965ea79SLois Curfman McInnes { 592dfbe8321SBarry Smith PetscErrorCode ierr; 59332077d6dSBarry Smith PetscTruth iascii,isbinary,isdraw,issocket; 5948965ea79SLois Curfman McInnes 595433994e6SBarry Smith PetscFunctionBegin; 5960f5bd95cSBarry Smith 59732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);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 60232077d6dSBarry Smith if (iascii || 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 { 607958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"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" 614dfbe8321SBarry Smith PetscErrorCode 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; 618dfbe8321SBarry Smith PetscErrorCode 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) { 637d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,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) { 644d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,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" 659dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op) 6608965ea79SLois Curfman McInnes { 66139ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 662dfbe8321SBarry Smith PetscErrorCode 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: 682b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 68312c028f9SKris Buschelman break; 68412c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 685273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 686273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 68712c028f9SKris Buschelman break; 68812c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 689273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 69012c028f9SKris Buschelman break; 69112c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 69229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 69377e54ba9SKris Buschelman case MAT_SYMMETRIC: 69477e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 6959a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 6969a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 6979a4540c5SBarry Smith case MAT_HERMITIAN: 6989a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 6999a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 7009a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 70177e54ba9SKris Buschelman break; 70212c028f9SKris Buschelman default: 70329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 7043a40ed3dSBarry Smith } 7053a40ed3dSBarry Smith PetscFunctionReturn(0); 7068965ea79SLois Curfman McInnes } 7078965ea79SLois Curfman McInnes 7084a2ae208SSatish Balay #undef __FUNCT__ 7094a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIDense" 710dfbe8321SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v) 7118965ea79SLois Curfman McInnes { 7123501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7136849ba73SBarry Smith PetscErrorCode ierr; 7146849ba73SBarry Smith int lrow,rstart = mat->rstart,rend = mat->rend; 7158965ea79SLois Curfman McInnes 7163a40ed3dSBarry Smith PetscFunctionBegin; 71729bbc08cSBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 7188965ea79SLois Curfman McInnes lrow = row - rstart; 719b3cc6726SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const int **)idx,(const PetscScalar **)v);CHKERRQ(ierr); 7203a40ed3dSBarry Smith PetscFunctionReturn(0); 7218965ea79SLois Curfman McInnes } 7228965ea79SLois Curfman McInnes 7234a2ae208SSatish Balay #undef __FUNCT__ 7244a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIDense" 725dfbe8321SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v) 7268965ea79SLois Curfman McInnes { 727dfbe8321SBarry Smith PetscErrorCode ierr; 728606d414cSSatish Balay 7293a40ed3dSBarry Smith PetscFunctionBegin; 730606d414cSSatish Balay if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 731606d414cSSatish Balay if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 7323a40ed3dSBarry Smith PetscFunctionReturn(0); 7338965ea79SLois Curfman McInnes } 7348965ea79SLois Curfman McInnes 7354a2ae208SSatish Balay #undef __FUNCT__ 7364a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 737dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7385b2fa520SLois Curfman McInnes { 7395b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7405b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 74187828ca2SBarry Smith PetscScalar *l,*r,x,*v; 742dfbe8321SBarry Smith PetscErrorCode ierr; 743dfbe8321SBarry Smith int i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 7445b2fa520SLois Curfman McInnes 7455b2fa520SLois Curfman McInnes PetscFunctionBegin; 74672d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7475b2fa520SLois Curfman McInnes if (ll) { 74872d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 74929bbc08cSBarry Smith if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 7501ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7515b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7525b2fa520SLois Curfman McInnes x = l[i]; 7535b2fa520SLois Curfman McInnes v = mat->v + i; 7545b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7555b2fa520SLois Curfman McInnes } 7561ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 757b0a32e0cSBarry Smith PetscLogFlops(n*m); 7585b2fa520SLois Curfman McInnes } 7595b2fa520SLois Curfman McInnes if (rr) { 76072d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 76129bbc08cSBarry Smith if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 7625b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7635b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7641ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 7655b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 7665b2fa520SLois Curfman McInnes x = r[i]; 7675b2fa520SLois Curfman McInnes v = mat->v + i*m; 7685b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 7695b2fa520SLois Curfman McInnes } 7701ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 771b0a32e0cSBarry Smith PetscLogFlops(n*m); 7725b2fa520SLois Curfman McInnes } 7735b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 7745b2fa520SLois Curfman McInnes } 7755b2fa520SLois Curfman McInnes 7764a2ae208SSatish Balay #undef __FUNCT__ 7774a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 778dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 779096963f5SLois Curfman McInnes { 7803501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7813501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 782dfbe8321SBarry Smith PetscErrorCode ierr; 783dfbe8321SBarry Smith int i,j; 784329f5518SBarry Smith PetscReal sum = 0.0; 78587828ca2SBarry Smith PetscScalar *v = mat->v; 7863501a2bdSLois Curfman McInnes 7873a40ed3dSBarry Smith PetscFunctionBegin; 7883501a2bdSLois Curfman McInnes if (mdn->size == 1) { 789064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 7903501a2bdSLois Curfman McInnes } else { 7913501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 792273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 793aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 794329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 7953501a2bdSLois Curfman McInnes #else 7963501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 7973501a2bdSLois Curfman McInnes #endif 7983501a2bdSLois Curfman McInnes } 799064f8208SBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 800064f8208SBarry Smith *nrm = sqrt(*nrm); 801b0a32e0cSBarry Smith PetscLogFlops(2*mdn->A->n*mdn->A->m); 8023a40ed3dSBarry Smith } else if (type == NORM_1) { 803329f5518SBarry Smith PetscReal *tmp,*tmp2; 804b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 805273d9f13SBarry Smith tmp2 = tmp + A->N; 806273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 807064f8208SBarry Smith *nrm = 0.0; 8083501a2bdSLois Curfman McInnes v = mat->v; 809273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 810273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 81167e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8123501a2bdSLois Curfman McInnes } 8133501a2bdSLois Curfman McInnes } 814d7d1e502SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 815273d9f13SBarry Smith for (j=0; j<A->N; j++) { 816064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8173501a2bdSLois Curfman McInnes } 818606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 819b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 8203a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 821329f5518SBarry Smith PetscReal ntemp; 8223501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 823064f8208SBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 8243a40ed3dSBarry Smith } else { 82529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8263501a2bdSLois Curfman McInnes } 8273501a2bdSLois Curfman McInnes } 8283a40ed3dSBarry Smith PetscFunctionReturn(0); 8293501a2bdSLois Curfman McInnes } 8303501a2bdSLois Curfman McInnes 8314a2ae208SSatish Balay #undef __FUNCT__ 8324a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 833dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout) 8343501a2bdSLois Curfman McInnes { 8353501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8363501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8373501a2bdSLois Curfman McInnes Mat B; 838273d9f13SBarry Smith int M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 8396849ba73SBarry Smith PetscErrorCode ierr; 8406849ba73SBarry Smith int j,i; 84187828ca2SBarry Smith PetscScalar *v; 8423501a2bdSLois Curfman McInnes 8433a40ed3dSBarry Smith PetscFunctionBegin; 8447c922b88SBarry Smith if (!matout && M != N) { 84529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 8467056b6fcSBarry Smith } 847878740d9SKris Buschelman ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B);CHKERRQ(ierr); 848878740d9SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 849878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);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" 872dfbe8321SBarry Smith PetscErrorCode MatScale_MPIDense(const 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; 8764ce68768SBarry Smith PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N; 87744cd7ae7SLois Curfman McInnes 8783a40ed3dSBarry Smith PetscFunctionBegin; 879268466fbSBarry Smith BLscal_(&nz,(PetscScalar*)alpha,a->v,&one); 880b0a32e0cSBarry Smith PetscLogFlops(nz); 8813a40ed3dSBarry Smith PetscFunctionReturn(0); 88244cd7ae7SLois Curfman McInnes } 88344cd7ae7SLois Curfman McInnes 8846849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 8858965ea79SLois Curfman McInnes 8864a2ae208SSatish Balay #undef __FUNCT__ 8874a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 888dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A) 889273d9f13SBarry Smith { 890dfbe8321SBarry Smith PetscErrorCode ierr; 891273d9f13SBarry Smith 892273d9f13SBarry Smith PetscFunctionBegin; 893273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 894273d9f13SBarry Smith PetscFunctionReturn(0); 895273d9f13SBarry Smith } 896273d9f13SBarry Smith 8978965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 89809dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 89909dc0095SBarry Smith MatGetRow_MPIDense, 90009dc0095SBarry Smith MatRestoreRow_MPIDense, 90109dc0095SBarry Smith MatMult_MPIDense, 90297304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 9037c922b88SBarry Smith MatMultTranspose_MPIDense, 9047c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 9058965ea79SLois Curfman McInnes 0, 90609dc0095SBarry Smith 0, 90709dc0095SBarry Smith 0, 90897304618SKris Buschelman /*10*/ 0, 90909dc0095SBarry Smith 0, 91009dc0095SBarry Smith 0, 91109dc0095SBarry Smith 0, 91209dc0095SBarry Smith MatTranspose_MPIDense, 91397304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 91497304618SKris Buschelman 0, 91509dc0095SBarry Smith MatGetDiagonal_MPIDense, 9165b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 91709dc0095SBarry Smith MatNorm_MPIDense, 91897304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 91909dc0095SBarry Smith MatAssemblyEnd_MPIDense, 92009dc0095SBarry Smith 0, 92109dc0095SBarry Smith MatSetOption_MPIDense, 92209dc0095SBarry Smith MatZeroEntries_MPIDense, 92397304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 92409dc0095SBarry Smith 0, 92509dc0095SBarry Smith 0, 92609dc0095SBarry Smith 0, 92709dc0095SBarry Smith 0, 92897304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 929273d9f13SBarry Smith 0, 93009dc0095SBarry Smith 0, 93109dc0095SBarry Smith MatGetArray_MPIDense, 93209dc0095SBarry Smith MatRestoreArray_MPIDense, 93397304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 93409dc0095SBarry Smith 0, 93509dc0095SBarry Smith 0, 93609dc0095SBarry Smith 0, 93709dc0095SBarry Smith 0, 93897304618SKris Buschelman /*40*/ 0, 9392ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 94009dc0095SBarry Smith 0, 94109dc0095SBarry Smith MatGetValues_MPIDense, 94209dc0095SBarry Smith 0, 94397304618SKris Buschelman /*45*/ 0, 94409dc0095SBarry Smith MatScale_MPIDense, 94509dc0095SBarry Smith 0, 94609dc0095SBarry Smith 0, 94709dc0095SBarry Smith 0, 948*521d7252SBarry Smith /*50*/ 0, 94909dc0095SBarry Smith 0, 95009dc0095SBarry Smith 0, 95109dc0095SBarry Smith 0, 95209dc0095SBarry Smith 0, 95397304618SKris Buschelman /*55*/ 0, 95409dc0095SBarry Smith 0, 95509dc0095SBarry Smith 0, 95609dc0095SBarry Smith 0, 95709dc0095SBarry Smith 0, 95897304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 959b9b97703SBarry Smith MatDestroy_MPIDense, 960b9b97703SBarry Smith MatView_MPIDense, 96197304618SKris Buschelman MatGetPetscMaps_Petsc, 96297304618SKris Buschelman 0, 96397304618SKris Buschelman /*65*/ 0, 96497304618SKris Buschelman 0, 96597304618SKris Buschelman 0, 96697304618SKris Buschelman 0, 96797304618SKris Buschelman 0, 96897304618SKris Buschelman /*70*/ 0, 96997304618SKris Buschelman 0, 97097304618SKris Buschelman 0, 97197304618SKris Buschelman 0, 97297304618SKris Buschelman 0, 97397304618SKris Buschelman /*75*/ 0, 97497304618SKris Buschelman 0, 97597304618SKris Buschelman 0, 97697304618SKris Buschelman 0, 97797304618SKris Buschelman 0, 97897304618SKris Buschelman /*80*/ 0, 97997304618SKris Buschelman 0, 98097304618SKris Buschelman 0, 98197304618SKris Buschelman 0, 982865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense, 983865e5f61SKris Buschelman 0, 984865e5f61SKris Buschelman 0, 985865e5f61SKris Buschelman 0, 986865e5f61SKris Buschelman 0, 987865e5f61SKris Buschelman 0, 988865e5f61SKris Buschelman /*90*/ 0, 989865e5f61SKris Buschelman 0, 990865e5f61SKris Buschelman 0, 991865e5f61SKris Buschelman 0, 992865e5f61SKris Buschelman 0, 993865e5f61SKris Buschelman /*95*/ 0, 994865e5f61SKris Buschelman 0, 995865e5f61SKris Buschelman 0, 996865e5f61SKris Buschelman 0}; 9978965ea79SLois Curfman McInnes 998273d9f13SBarry Smith EXTERN_C_BEGIN 9994a2ae208SSatish Balay #undef __FUNCT__ 1000a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 1001dfbe8321SBarry Smith PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1002a23d5eceSKris Buschelman { 1003a23d5eceSKris Buschelman Mat_MPIDense *a; 1004dfbe8321SBarry Smith PetscErrorCode ierr; 1005a23d5eceSKris Buschelman 1006a23d5eceSKris Buschelman PetscFunctionBegin; 1007a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1008a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1009a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1010a23d5eceSKris Buschelman 1011a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 10125c5985e7SKris Buschelman ierr = MatCreate(PETSC_COMM_SELF,mat->m,mat->N,mat->m,mat->N,&a->A);CHKERRQ(ierr); 10135c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 10145c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 1015a23d5eceSKris Buschelman PetscLogObjectParent(mat,a->A); 1016a23d5eceSKris Buschelman PetscFunctionReturn(0); 1017a23d5eceSKris Buschelman } 1018a23d5eceSKris Buschelman EXTERN_C_END 1019a23d5eceSKris Buschelman 10200bad9183SKris Buschelman /*MC 1021fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 10220bad9183SKris Buschelman 10230bad9183SKris Buschelman Options Database Keys: 10240bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 10250bad9183SKris Buschelman 10260bad9183SKris Buschelman Level: beginner 10270bad9183SKris Buschelman 10280bad9183SKris Buschelman .seealso: MatCreateMPIDense 10290bad9183SKris Buschelman M*/ 10300bad9183SKris Buschelman 1031a23d5eceSKris Buschelman EXTERN_C_BEGIN 1032a23d5eceSKris Buschelman #undef __FUNCT__ 10334a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1034dfbe8321SBarry Smith PetscErrorCode MatCreate_MPIDense(Mat mat) 1035273d9f13SBarry Smith { 1036273d9f13SBarry Smith Mat_MPIDense *a; 1037dfbe8321SBarry Smith PetscErrorCode ierr; 1038dfbe8321SBarry Smith int i; 1039273d9f13SBarry Smith 1040273d9f13SBarry Smith PetscFunctionBegin; 1041b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1042b0a32e0cSBarry Smith mat->data = (void*)a; 1043273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1044273d9f13SBarry Smith mat->factor = 0; 1045273d9f13SBarry Smith mat->mapping = 0; 1046273d9f13SBarry Smith 1047273d9f13SBarry Smith a->factor = 0; 1048273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1049273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1050273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1051273d9f13SBarry Smith 1052273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1053273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1054273d9f13SBarry Smith a->nvec = mat->n; 1055273d9f13SBarry Smith 1056273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 1057273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 10588a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 10598a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1060273d9f13SBarry Smith 1061273d9f13SBarry Smith /* build local table of row and column ownerships */ 1062b0a32e0cSBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1063273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 1064b0a32e0cSBarry Smith PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1065273d9f13SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1066273d9f13SBarry Smith a->rowners[0] = 0; 1067273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1068273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1069273d9f13SBarry Smith } 1070273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1071273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 1072273d9f13SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1073273d9f13SBarry Smith a->cowners[0] = 0; 1074273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1075273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1076273d9f13SBarry Smith } 1077273d9f13SBarry Smith 1078273d9f13SBarry Smith /* build cache for off array entries formed */ 1079273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1080273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1081273d9f13SBarry Smith 1082273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1083273d9f13SBarry Smith a->lvec = 0; 1084273d9f13SBarry Smith a->Mvctx = 0; 1085273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1086273d9f13SBarry Smith 1087273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1088273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1089273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1090a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1091a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1092a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1093273d9f13SBarry Smith PetscFunctionReturn(0); 1094273d9f13SBarry Smith } 1095273d9f13SBarry Smith EXTERN_C_END 1096273d9f13SBarry Smith 1097209238afSKris Buschelman /*MC 1098002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1099209238afSKris Buschelman 1100209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1101209238afSKris Buschelman and MATMPIDENSE otherwise. 1102209238afSKris Buschelman 1103209238afSKris Buschelman Options Database Keys: 1104209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1105209238afSKris Buschelman 1106209238afSKris Buschelman Level: beginner 1107209238afSKris Buschelman 1108209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1109209238afSKris Buschelman M*/ 1110209238afSKris Buschelman 1111209238afSKris Buschelman EXTERN_C_BEGIN 1112209238afSKris Buschelman #undef __FUNCT__ 1113209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1114dfbe8321SBarry Smith PetscErrorCode MatCreate_Dense(Mat A) 1115dfbe8321SBarry Smith { 11166849ba73SBarry Smith PetscErrorCode ierr; 11176849ba73SBarry Smith int size; 1118209238afSKris Buschelman 1119209238afSKris Buschelman PetscFunctionBegin; 1120209238afSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1121209238afSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1122209238afSKris Buschelman if (size == 1) { 1123209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1124209238afSKris Buschelman } else { 1125209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1126209238afSKris Buschelman } 1127209238afSKris Buschelman PetscFunctionReturn(0); 1128209238afSKris Buschelman } 1129209238afSKris Buschelman EXTERN_C_END 1130209238afSKris Buschelman 11314a2ae208SSatish Balay #undef __FUNCT__ 11324a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1133273d9f13SBarry Smith /*@C 1134273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1135273d9f13SBarry Smith 1136273d9f13SBarry Smith Not collective 1137273d9f13SBarry Smith 1138273d9f13SBarry Smith Input Parameters: 1139273d9f13SBarry Smith . A - the matrix 1140273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1141273d9f13SBarry Smith to control all matrix memory allocation. 1142273d9f13SBarry Smith 1143273d9f13SBarry Smith Notes: 1144273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1145273d9f13SBarry Smith storage by columns. 1146273d9f13SBarry Smith 1147273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1148273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1149273d9f13SBarry Smith set data=PETSC_NULL. 1150273d9f13SBarry Smith 1151273d9f13SBarry Smith Level: intermediate 1152273d9f13SBarry Smith 1153273d9f13SBarry Smith .keywords: matrix,dense, parallel 1154273d9f13SBarry Smith 1155273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1156273d9f13SBarry Smith @*/ 1157dfbe8321SBarry Smith PetscErrorCode MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1158273d9f13SBarry Smith { 1159dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar *); 1160273d9f13SBarry Smith 1161273d9f13SBarry Smith PetscFunctionBegin; 1162565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1163a23d5eceSKris Buschelman if (f) { 1164a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1165a23d5eceSKris Buschelman } 1166273d9f13SBarry Smith PetscFunctionReturn(0); 1167273d9f13SBarry Smith } 1168273d9f13SBarry Smith 11694a2ae208SSatish Balay #undef __FUNCT__ 11704a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 11718965ea79SLois Curfman McInnes /*@C 117239ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 11738965ea79SLois Curfman McInnes 1174db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1175db81eaa0SLois Curfman McInnes 11768965ea79SLois Curfman McInnes Input Parameters: 1177db81eaa0SLois Curfman McInnes + comm - MPI communicator 11788965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1179db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 11808965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1181db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 11827f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1183dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 11848965ea79SLois Curfman McInnes 11858965ea79SLois Curfman McInnes Output Parameter: 1186477f1c0bSLois Curfman McInnes . A - the matrix 11878965ea79SLois Curfman McInnes 1188b259b22eSLois Curfman McInnes Notes: 118939ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 119039ddd567SLois Curfman McInnes storage by columns. 11918965ea79SLois Curfman McInnes 119218f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 119318f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 11947f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 119518f449edSLois Curfman McInnes 11968965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 11978965ea79SLois Curfman McInnes (possibly both). 11988965ea79SLois Curfman McInnes 1199027ccd11SLois Curfman McInnes Level: intermediate 1200027ccd11SLois Curfman McInnes 120139ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 12028965ea79SLois Curfman McInnes 120339ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 12048965ea79SLois Curfman McInnes @*/ 1205dfbe8321SBarry Smith PetscErrorCode MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A) 12068965ea79SLois Curfman McInnes { 12076849ba73SBarry Smith PetscErrorCode ierr; 12086849ba73SBarry Smith int size; 12098965ea79SLois Curfman McInnes 12103a40ed3dSBarry Smith PetscFunctionBegin; 1211273d9f13SBarry Smith ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr); 1212273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1213273d9f13SBarry Smith if (size > 1) { 1214273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1215273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1216273d9f13SBarry Smith } else { 1217273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1218273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 12198c469469SLois Curfman McInnes } 12203a40ed3dSBarry Smith PetscFunctionReturn(0); 12218965ea79SLois Curfman McInnes } 12228965ea79SLois Curfman McInnes 12234a2ae208SSatish Balay #undef __FUNCT__ 12244a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 12256849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 12268965ea79SLois Curfman McInnes { 12278965ea79SLois Curfman McInnes Mat mat; 12283501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1229dfbe8321SBarry Smith PetscErrorCode ierr; 12308965ea79SLois Curfman McInnes 12313a40ed3dSBarry Smith PetscFunctionBegin; 12328965ea79SLois Curfman McInnes *newmat = 0; 1233273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr); 1234be5d1d56SKris Buschelman ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr); 1235b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1236b0a32e0cSBarry Smith mat->data = (void*)a; 1237549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 12383501a2bdSLois Curfman McInnes mat->factor = A->factor; 1239c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1240273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 12418965ea79SLois Curfman McInnes 12428965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 12438965ea79SLois Curfman McInnes a->rend = oldmat->rend; 12448965ea79SLois Curfman McInnes a->size = oldmat->size; 12458965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1246e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1247b9b97703SBarry Smith a->nvec = oldmat->nvec; 12483782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1249b0a32e0cSBarry Smith ierr = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1250b0a32e0cSBarry Smith PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1251549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 12528798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 12538965ea79SLois Curfman McInnes 1254329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 12555609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 1256b0a32e0cSBarry Smith PetscLogObjectParent(mat,a->A); 12578965ea79SLois Curfman McInnes *newmat = mat; 12583a40ed3dSBarry Smith PetscFunctionReturn(0); 12598965ea79SLois Curfman McInnes } 12608965ea79SLois Curfman McInnes 1261e090d566SSatish Balay #include "petscsys.h" 12628965ea79SLois Curfman McInnes 12634a2ae208SSatish Balay #undef __FUNCT__ 12644a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 1265dfbe8321SBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,const MatType type,Mat *newmat) 126690ace30eSBarry Smith { 12676849ba73SBarry Smith PetscErrorCode ierr; 126832dcc486SBarry Smith PetscMPIInt rank,size; 126932dcc486SBarry Smith int *rowners,i,m,nz,j; 127087828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 127190ace30eSBarry Smith MPI_Status status; 127290ace30eSBarry Smith 12733a40ed3dSBarry Smith PetscFunctionBegin; 1274d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1275d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 127690ace30eSBarry Smith 127790ace30eSBarry Smith /* determine ownership of all rows */ 127890ace30eSBarry Smith m = M/size + ((M % size) > rank); 1279b0a32e0cSBarry Smith ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1280ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 128190ace30eSBarry Smith rowners[0] = 0; 128290ace30eSBarry Smith for (i=2; i<=size; i++) { 128390ace30eSBarry Smith rowners[i] += rowners[i-1]; 128490ace30eSBarry Smith } 128590ace30eSBarry Smith 1286878740d9SKris Buschelman ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 1287be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1288878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 128990ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 129090ace30eSBarry Smith 129190ace30eSBarry Smith if (!rank) { 129287828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 129390ace30eSBarry Smith 129490ace30eSBarry Smith /* read in my part of the matrix numerical values */ 12950752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 129690ace30eSBarry Smith 129790ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 129890ace30eSBarry Smith vals_ptr = vals; 129990ace30eSBarry Smith for (i=0; i<m; i++) { 130090ace30eSBarry Smith for (j=0; j<N; j++) { 130190ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 130290ace30eSBarry Smith } 130390ace30eSBarry Smith } 130490ace30eSBarry Smith 130590ace30eSBarry Smith /* read in other processors and ship out */ 130690ace30eSBarry Smith for (i=1; i<size; i++) { 130790ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 13080752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1309ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 131090ace30eSBarry Smith } 13113a40ed3dSBarry Smith } else { 131290ace30eSBarry Smith /* receive numeric values */ 131387828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 131490ace30eSBarry Smith 131590ace30eSBarry Smith /* receive message of values*/ 1316ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 131790ace30eSBarry Smith 131890ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 131990ace30eSBarry Smith vals_ptr = vals; 132090ace30eSBarry Smith for (i=0; i<m; i++) { 132190ace30eSBarry Smith for (j=0; j<N; j++) { 132290ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 132390ace30eSBarry Smith } 132490ace30eSBarry Smith } 132590ace30eSBarry Smith } 1326606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1327606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 13286d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13296d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13303a40ed3dSBarry Smith PetscFunctionReturn(0); 133190ace30eSBarry Smith } 133290ace30eSBarry Smith 13334a2ae208SSatish Balay #undef __FUNCT__ 13344a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1335dfbe8321SBarry Smith PetscErrorCode MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat) 13368965ea79SLois Curfman McInnes { 13378965ea79SLois Curfman McInnes Mat A; 133887828ca2SBarry Smith PetscScalar *vals,*svals; 133919bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 13408965ea79SLois Curfman McInnes MPI_Status status; 134132dcc486SBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag; 134232dcc486SBarry Smith int header[4],*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 13438965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols; 13446849ba73SBarry Smith int i,nz,j,rstart,rend,fd; 13456849ba73SBarry Smith PetscErrorCode ierr; 13468965ea79SLois Curfman McInnes 13473a40ed3dSBarry Smith PetscFunctionBegin; 1348d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1349d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13508965ea79SLois Curfman McInnes if (!rank) { 1351b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 13520752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1353552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 13548965ea79SLois Curfman McInnes } 13558965ea79SLois Curfman McInnes 1356ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 135790ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 135890ace30eSBarry Smith 135990ace30eSBarry Smith /* 136090ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 136190ace30eSBarry Smith */ 136290ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1363be5d1d56SKris Buschelman ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 13643a40ed3dSBarry Smith PetscFunctionReturn(0); 136590ace30eSBarry Smith } 136690ace30eSBarry Smith 13678965ea79SLois Curfman McInnes /* determine ownership of all rows */ 13688965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 1369b0a32e0cSBarry Smith ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1370ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 13718965ea79SLois Curfman McInnes rowners[0] = 0; 13728965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 13738965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 13748965ea79SLois Curfman McInnes } 13758965ea79SLois Curfman McInnes rstart = rowners[rank]; 13768965ea79SLois Curfman McInnes rend = rowners[rank+1]; 13778965ea79SLois Curfman McInnes 13788965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 137984d528b1SBarry Smith ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(int),&ourlens);CHKERRQ(ierr); 13808965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 13818965ea79SLois Curfman McInnes if (!rank) { 1382b0a32e0cSBarry Smith ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr); 13830752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 1384b0a32e0cSBarry Smith ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr); 13858965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 1386ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1387606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1388ca161407SBarry Smith } else { 1389ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 13908965ea79SLois Curfman McInnes } 13918965ea79SLois Curfman McInnes 13928965ea79SLois Curfman McInnes if (!rank) { 13938965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 1394b0a32e0cSBarry Smith ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr); 1395549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 13968965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 13978965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 13988965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 13998965ea79SLois Curfman McInnes } 14008965ea79SLois Curfman McInnes } 1401606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 14028965ea79SLois Curfman McInnes 14038965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 14048965ea79SLois Curfman McInnes maxnz = 0; 14058965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14060452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 14078965ea79SLois Curfman McInnes } 1408b0a32e0cSBarry Smith ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr); 14098965ea79SLois Curfman McInnes 14108965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 14118965ea79SLois Curfman McInnes nz = procsnz[0]; 1412b0a32e0cSBarry Smith ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr); 14130752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 14148965ea79SLois Curfman McInnes 14158965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 14168965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14178965ea79SLois Curfman McInnes nz = procsnz[i]; 14180752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1419ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 14208965ea79SLois Curfman McInnes } 1421606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 14223a40ed3dSBarry Smith } else { 14238965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 14248965ea79SLois Curfman McInnes nz = 0; 14258965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14268965ea79SLois Curfman McInnes nz += ourlens[i]; 14278965ea79SLois Curfman McInnes } 142884d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(int),&mycols);CHKERRQ(ierr); 14298965ea79SLois Curfman McInnes 14308965ea79SLois Curfman McInnes /* receive message of column indices*/ 1431ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1432ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 143329bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14348965ea79SLois Curfman McInnes } 14358965ea79SLois Curfman McInnes 14368965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1437549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 14388965ea79SLois Curfman McInnes jj = 0; 14398965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14408965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 14418965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 14428965ea79SLois Curfman McInnes jj++; 14438965ea79SLois Curfman McInnes } 14448965ea79SLois Curfman McInnes } 14458965ea79SLois Curfman McInnes 14468965ea79SLois Curfman McInnes /* create our matrix */ 14478965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14488965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 14498965ea79SLois Curfman McInnes } 1450878740d9SKris Buschelman ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 1451be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1452878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 14538965ea79SLois Curfman McInnes A = *newmat; 14548965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14558965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 14568965ea79SLois Curfman McInnes } 14578965ea79SLois Curfman McInnes 14588965ea79SLois Curfman McInnes if (!rank) { 145987828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14608965ea79SLois Curfman McInnes 14618965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 14628965ea79SLois Curfman McInnes nz = procsnz[0]; 14630752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 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 14768965ea79SLois Curfman McInnes /* read in other processors and ship out */ 14778965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14788965ea79SLois Curfman McInnes nz = procsnz[i]; 14790752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1480ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 14818965ea79SLois Curfman McInnes } 1482606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 14833a40ed3dSBarry Smith } else { 14848965ea79SLois Curfman McInnes /* receive numeric values */ 148584d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14868965ea79SLois Curfman McInnes 14878965ea79SLois Curfman McInnes /* receive message of values*/ 1488ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1489ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 149029bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14918965ea79SLois Curfman McInnes 14928965ea79SLois Curfman McInnes /* insert into matrix */ 14938965ea79SLois Curfman McInnes jj = rstart; 14948965ea79SLois Curfman McInnes smycols = mycols; 14958965ea79SLois Curfman McInnes svals = vals; 14968965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14978965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14988965ea79SLois Curfman McInnes smycols += ourlens[i]; 14998965ea79SLois Curfman McInnes svals += ourlens[i]; 15008965ea79SLois Curfman McInnes jj++; 15018965ea79SLois Curfman McInnes } 15028965ea79SLois Curfman McInnes } 1503606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1504606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1505606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1506606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 15078965ea79SLois Curfman McInnes 15086d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15096d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15103a40ed3dSBarry Smith PetscFunctionReturn(0); 15118965ea79SLois Curfman McInnes } 151290ace30eSBarry Smith 151390ace30eSBarry Smith 151490ace30eSBarry Smith 151590ace30eSBarry Smith 151690ace30eSBarry Smith 1517