173f4d377SMatthew Knepley /*$Id: mpidense.c,v 1.159 2001/08/10 03:30:41 bsmith Exp $*/ 28965ea79SLois Curfman McInnes 3ed3cc1f0SBarry Smith /* 4ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 5ed3cc1f0SBarry Smith */ 6ed3cc1f0SBarry Smith 770f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 8f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 98965ea79SLois Curfman McInnes 100de54da6SSatish Balay EXTERN_C_BEGIN 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 130de54da6SSatish Balay int MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 140de54da6SSatish Balay { 150de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 16cfce73b9SSatish Balay int m = A->m,rstart = mdn->rstart,ierr; 1787828ca2SBarry Smith PetscScalar *array; 180de54da6SSatish Balay MPI_Comm comm; 190de54da6SSatish Balay 200de54da6SSatish Balay PetscFunctionBegin; 21273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 220de54da6SSatish Balay 230de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 240de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 250de54da6SSatish Balay 260de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 270de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 280de54da6SSatish Balay ierr = MatCreateSeqDense(comm,m,m,array+m*rstart,B);CHKERRQ(ierr); 290de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 300de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 310de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 320de54da6SSatish Balay 330de54da6SSatish Balay *iscopy = PETSC_TRUE; 340de54da6SSatish Balay PetscFunctionReturn(0); 350de54da6SSatish Balay } 360de54da6SSatish Balay EXTERN_C_END 370de54da6SSatish Balay 384a2ae208SSatish Balay #undef __FUNCT__ 394a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 40f15d580aSBarry Smith int 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; 4339b7565bSBarry Smith int ierr,i,j,rstart = A->rstart,rend = A->rend,row; 44273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 458965ea79SLois Curfman McInnes 463a40ed3dSBarry Smith PetscFunctionBegin; 478965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 485ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 49273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 508965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 518965ea79SLois Curfman McInnes row = idxm[i] - rstart; 5239b7565bSBarry Smith if (roworiented) { 5339b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 543a40ed3dSBarry Smith } else { 558965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 565ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 57273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 5839b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 5939b7565bSBarry Smith } 608965ea79SLois Curfman McInnes } 613a40ed3dSBarry Smith } else { 623782ba37SSatish Balay if (!A->donotstash) { 6339b7565bSBarry Smith if (roworiented) { 648798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 65d36fbae8SSatish Balay } else { 668798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 6739b7565bSBarry Smith } 68b49de8d1SLois Curfman McInnes } 69b49de8d1SLois Curfman McInnes } 703782ba37SSatish Balay } 713a40ed3dSBarry Smith PetscFunctionReturn(0); 72b49de8d1SLois Curfman McInnes } 73b49de8d1SLois Curfman McInnes 744a2ae208SSatish Balay #undef __FUNCT__ 754a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 76f15d580aSBarry Smith int MatGetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],PetscScalar v[]) 77b49de8d1SLois Curfman McInnes { 78b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 79b49de8d1SLois Curfman McInnes int ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row; 80b49de8d1SLois Curfman McInnes 813a40ed3dSBarry Smith PetscFunctionBegin; 82b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 8329bbc08cSBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 84273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 85b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 86b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 87b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 8829bbc08cSBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 89273d9f13SBarry Smith if (idxn[j] >= mat->N) { 9029bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 91a8c6a408SBarry Smith } 92b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 93b49de8d1SLois Curfman McInnes } 94a8c6a408SBarry Smith } else { 9529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 968965ea79SLois Curfman McInnes } 978965ea79SLois Curfman McInnes } 983a40ed3dSBarry Smith PetscFunctionReturn(0); 998965ea79SLois Curfman McInnes } 1008965ea79SLois Curfman McInnes 1014a2ae208SSatish Balay #undef __FUNCT__ 1024a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 1034e7234bfSBarry Smith int MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 104ff14e315SSatish Balay { 105ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 106ff14e315SSatish Balay int ierr; 107ff14e315SSatish Balay 1083a40ed3dSBarry Smith PetscFunctionBegin; 109ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1103a40ed3dSBarry Smith PetscFunctionReturn(0); 111ff14e315SSatish Balay } 112ff14e315SSatish Balay 1134a2ae208SSatish Balay #undef __FUNCT__ 1144a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 115ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 116ca3fa75bSLois Curfman McInnes { 117ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 118ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 119cfce73b9SSatish Balay int i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 12087828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 121ca3fa75bSLois Curfman McInnes Mat newmat; 122ca3fa75bSLois Curfman McInnes 123ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 124ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 125ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 126b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 127b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 128ca3fa75bSLois Curfman McInnes 129ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1307eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1317eba5e9cSLois Curfman McInnes original matrix! */ 132ca3fa75bSLois Curfman McInnes 133ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1347eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 135ca3fa75bSLois Curfman McInnes 136ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 137ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 13829bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 1397eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 140ca3fa75bSLois Curfman McInnes newmat = *B; 141ca3fa75bSLois Curfman McInnes } else { 142ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 14332828cfdSBarry Smith ierr = MatCreateMPIDense(A->comm,nrows,cs,PETSC_DECIDE,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 144ca3fa75bSLois Curfman McInnes } 145ca3fa75bSLois Curfman McInnes 146ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 147ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 148ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 149ca3fa75bSLois Curfman McInnes 150ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 151ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 152ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 1537eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 154ca3fa75bSLois Curfman McInnes } 155ca3fa75bSLois Curfman McInnes } 156ca3fa75bSLois Curfman McInnes 157ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 158ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 159ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 160ca3fa75bSLois Curfman McInnes 161ca3fa75bSLois Curfman McInnes /* Free work space */ 162ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 163ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 164ca3fa75bSLois Curfman McInnes *B = newmat; 165ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 166ca3fa75bSLois Curfman McInnes } 167ca3fa75bSLois Curfman McInnes 1684a2ae208SSatish Balay #undef __FUNCT__ 1694a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 1704e7234bfSBarry Smith int MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 171ff14e315SSatish Balay { 1723a40ed3dSBarry Smith PetscFunctionBegin; 1733a40ed3dSBarry Smith PetscFunctionReturn(0); 174ff14e315SSatish Balay } 175ff14e315SSatish Balay 1764a2ae208SSatish Balay #undef __FUNCT__ 1774a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 1788f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 1798965ea79SLois Curfman McInnes { 18039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 1818965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 182d36fbae8SSatish Balay int ierr,nstash,reallocs; 1838965ea79SLois Curfman McInnes InsertMode addv; 1848965ea79SLois Curfman McInnes 1853a40ed3dSBarry Smith PetscFunctionBegin; 1868965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 187ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1887056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 18929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 1908965ea79SLois Curfman McInnes } 191e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1928965ea79SLois Curfman McInnes 1938798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 1948798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 195b0a32e0cSBarry Smith PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs); 1963a40ed3dSBarry Smith PetscFunctionReturn(0); 1978965ea79SLois Curfman McInnes } 1988965ea79SLois Curfman McInnes 1994a2ae208SSatish Balay #undef __FUNCT__ 2004a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 2018f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2028965ea79SLois Curfman McInnes { 20339ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2047ef1d9bdSSatish Balay int i,n,ierr,*row,*col,flg,j,rstart,ncols; 20587828ca2SBarry Smith PetscScalar *val; 206e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2078965ea79SLois Curfman McInnes 2083a40ed3dSBarry Smith PetscFunctionBegin; 2098965ea79SLois Curfman McInnes /* wait on receives */ 2107ef1d9bdSSatish Balay while (1) { 2118798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2127ef1d9bdSSatish Balay if (!flg) break; 2138965ea79SLois Curfman McInnes 2147ef1d9bdSSatish Balay for (i=0; i<n;) { 2157ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2167ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2177ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2187ef1d9bdSSatish Balay else ncols = n-i; 2197ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2207ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2217ef1d9bdSSatish Balay i = j; 2228965ea79SLois Curfman McInnes } 2237ef1d9bdSSatish Balay } 2248798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2258965ea79SLois Curfman McInnes 22639ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 22739ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2288965ea79SLois Curfman McInnes 2296d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 23039ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2318965ea79SLois Curfman McInnes } 2323a40ed3dSBarry Smith PetscFunctionReturn(0); 2338965ea79SLois Curfman McInnes } 2348965ea79SLois Curfman McInnes 2354a2ae208SSatish Balay #undef __FUNCT__ 2364a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 2378f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A) 2388965ea79SLois Curfman McInnes { 2393a40ed3dSBarry Smith int ierr; 24039ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2413a40ed3dSBarry Smith 2423a40ed3dSBarry Smith PetscFunctionBegin; 2433a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2443a40ed3dSBarry Smith PetscFunctionReturn(0); 2458965ea79SLois Curfman McInnes } 2468965ea79SLois Curfman McInnes 2474a2ae208SSatish Balay #undef __FUNCT__ 2484a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIDense" 2498f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs) 2504e220ebcSLois Curfman McInnes { 2513a40ed3dSBarry Smith PetscFunctionBegin; 2524e220ebcSLois Curfman McInnes *bs = 1; 2533a40ed3dSBarry Smith PetscFunctionReturn(0); 2544e220ebcSLois Curfman McInnes } 2554e220ebcSLois Curfman McInnes 2568965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2578965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2588965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2593501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2608965ea79SLois Curfman McInnes routine. 2618965ea79SLois Curfman McInnes */ 2624a2ae208SSatish Balay #undef __FUNCT__ 2634a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 264268466fbSBarry Smith int MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag) 2658965ea79SLois Curfman McInnes { 26639ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2678965ea79SLois Curfman McInnes int i,ierr,N,*rows,*owners = l->rowners,size = l->size; 268c1dc657dSBarry Smith int *nprocs,j,idx,nsends; 2698965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2708965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2718965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2728965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2738965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2748965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2758965ea79SLois Curfman McInnes IS istmp; 27635d8aa7fSBarry Smith PetscTruth found; 2778965ea79SLois Curfman McInnes 2783a40ed3dSBarry Smith PetscFunctionBegin; 279b9b97703SBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 2808965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 2818965ea79SLois Curfman McInnes 2828965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 283b0a32e0cSBarry Smith ierr = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr); 284549d3d68SSatish Balay ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 285b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/ 2868965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 2878965ea79SLois Curfman McInnes idx = rows[i]; 28835d8aa7fSBarry Smith found = PETSC_FALSE; 2898965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 2908965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 291c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 2928965ea79SLois Curfman McInnes } 2938965ea79SLois Curfman McInnes } 29429bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 2958965ea79SLois Curfman McInnes } 296c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 2978965ea79SLois Curfman McInnes 2988965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 299c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3008965ea79SLois Curfman McInnes 3018965ea79SLois Curfman McInnes /* post receives: */ 302b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr); 303b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3048965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 305ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3068965ea79SLois Curfman McInnes } 3078965ea79SLois Curfman McInnes 3088965ea79SLois Curfman McInnes /* do sends: 3098965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3108965ea79SLois Curfman McInnes the ith processor 3118965ea79SLois Curfman McInnes */ 312b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr); 313b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 314b0a32e0cSBarry Smith ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr); 3158965ea79SLois Curfman McInnes starts[0] = 0; 316c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3178965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3188965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3198965ea79SLois Curfman McInnes } 3208965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3218965ea79SLois Curfman McInnes 3228965ea79SLois Curfman McInnes starts[0] = 0; 323c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3248965ea79SLois Curfman McInnes count = 0; 3258965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 326c1dc657dSBarry Smith if (nprocs[2*i+1]) { 327c1dc657dSBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3288965ea79SLois Curfman McInnes } 3298965ea79SLois Curfman McInnes } 330606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3318965ea79SLois Curfman McInnes 3328965ea79SLois Curfman McInnes base = owners[rank]; 3338965ea79SLois Curfman McInnes 3348965ea79SLois Curfman McInnes /* wait on receives */ 335b0a32e0cSBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr); 3368965ea79SLois Curfman McInnes source = lens + nrecvs; 3378965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3388965ea79SLois Curfman McInnes while (count) { 339ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3408965ea79SLois Curfman McInnes /* unpack receives into our local space */ 341ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 3428965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3438965ea79SLois Curfman McInnes lens[imdex] = n; 3448965ea79SLois Curfman McInnes slen += n; 3458965ea79SLois Curfman McInnes count--; 3468965ea79SLois Curfman McInnes } 347606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3488965ea79SLois Curfman McInnes 3498965ea79SLois Curfman McInnes /* move the data into the send scatter */ 350b0a32e0cSBarry Smith ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr); 3518965ea79SLois Curfman McInnes count = 0; 3528965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 3538965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3548965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 3558965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3568965ea79SLois Curfman McInnes } 3578965ea79SLois Curfman McInnes } 358606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 359606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 360606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 361606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3628965ea79SLois Curfman McInnes 3638965ea79SLois Curfman McInnes /* actually zap the local rows */ 364029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 365b0a32e0cSBarry Smith PetscLogObjectParent(A,istmp); 366606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3678965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 3688965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 3698965ea79SLois Curfman McInnes 3708965ea79SLois Curfman McInnes /* wait on sends */ 3718965ea79SLois Curfman McInnes if (nsends) { 372b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 373ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 374606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 3758965ea79SLois Curfman McInnes } 376606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 377606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 3788965ea79SLois Curfman McInnes 3793a40ed3dSBarry Smith PetscFunctionReturn(0); 3808965ea79SLois Curfman McInnes } 3818965ea79SLois Curfman McInnes 3824a2ae208SSatish Balay #undef __FUNCT__ 3834a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 3848f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3858965ea79SLois Curfman McInnes { 38639ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 3878965ea79SLois Curfman McInnes int ierr; 388c456f294SBarry Smith 3893a40ed3dSBarry Smith PetscFunctionBegin; 39043a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39143a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39244cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 3933a40ed3dSBarry Smith PetscFunctionReturn(0); 3948965ea79SLois Curfman McInnes } 3958965ea79SLois Curfman McInnes 3964a2ae208SSatish Balay #undef __FUNCT__ 3974a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 3988f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3998965ea79SLois Curfman McInnes { 40039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4018965ea79SLois Curfman McInnes int ierr; 402c456f294SBarry Smith 4033a40ed3dSBarry Smith PetscFunctionBegin; 40443a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 40543a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 40644cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4073a40ed3dSBarry Smith PetscFunctionReturn(0); 4088965ea79SLois Curfman McInnes } 4098965ea79SLois Curfman McInnes 4104a2ae208SSatish Balay #undef __FUNCT__ 4114a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 4127c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 413096963f5SLois Curfman McInnes { 414096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 415096963f5SLois Curfman McInnes int ierr; 41687828ca2SBarry Smith PetscScalar zero = 0.0; 417096963f5SLois Curfman McInnes 4183a40ed3dSBarry Smith PetscFunctionBegin; 4193501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 4207c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 421537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 422537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4233a40ed3dSBarry Smith PetscFunctionReturn(0); 424096963f5SLois Curfman McInnes } 425096963f5SLois Curfman McInnes 4264a2ae208SSatish Balay #undef __FUNCT__ 4274a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 4287c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 429096963f5SLois Curfman McInnes { 430096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 431096963f5SLois Curfman McInnes int ierr; 432096963f5SLois Curfman McInnes 4333a40ed3dSBarry Smith PetscFunctionBegin; 4343501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4357c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 436537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 437537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4383a40ed3dSBarry Smith PetscFunctionReturn(0); 439096963f5SLois Curfman McInnes } 440096963f5SLois Curfman McInnes 4414a2ae208SSatish Balay #undef __FUNCT__ 4424a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 4438f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 4448965ea79SLois Curfman McInnes { 44539ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 446096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 447273d9f13SBarry Smith int ierr,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); 452096963f5SLois Curfman McInnes 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 } 4609a8c540fSSatish Balay 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" 466e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 4678965ea79SLois Curfman McInnes { 4683501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4698965ea79SLois Curfman McInnes int ierr; 470ed3cc1f0SBarry Smith 4713a40ed3dSBarry Smith PetscFunctionBegin; 47294d884c6SBarry Smith 473aa482453SBarry Smith #if defined(PETSC_USE_LOG) 474b0a32e0cSBarry 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); 4883a40ed3dSBarry Smith PetscFunctionReturn(0); 4898965ea79SLois Curfman McInnes } 49039ddd567SLois Curfman McInnes 4914a2ae208SSatish Balay #undef __FUNCT__ 4924a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 493b0a32e0cSBarry Smith static int MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 4948965ea79SLois Curfman McInnes { 49539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 4968965ea79SLois Curfman McInnes int ierr; 4977056b6fcSBarry Smith 4983a40ed3dSBarry Smith PetscFunctionBegin; 49939ddd567SLois Curfman McInnes if (mdn->size == 1) { 50039ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5018965ea79SLois Curfman McInnes } 50229bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5033a40ed3dSBarry Smith PetscFunctionReturn(0); 5048965ea79SLois Curfman McInnes } 5058965ea79SLois Curfman McInnes 5064a2ae208SSatish Balay #undef __FUNCT__ 5074a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 508b0a32e0cSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5098965ea79SLois Curfman McInnes { 51039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 511fb9695e5SSatish Balay int ierr,size = mdn->size,rank = mdn->rank; 512b0a32e0cSBarry Smith PetscViewerType vtype; 513f1af5d2fSBarry Smith PetscTruth isascii,isdraw; 514b0a32e0cSBarry Smith PetscViewer sviewer; 515f3ef73ceSBarry Smith PetscViewerFormat format; 5168965ea79SLois Curfman McInnes 5173a40ed3dSBarry Smith PetscFunctionBegin; 518b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 519fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 520f1af5d2fSBarry Smith if (isascii) { 521b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 522b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 523456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 5244e220ebcSLois Curfman McInnes MatInfo info; 525888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 526b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m, 5276831982aSBarry Smith (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr); 528b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5293501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5303a40ed3dSBarry Smith PetscFunctionReturn(0); 531fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 5323a40ed3dSBarry Smith PetscFunctionReturn(0); 5338965ea79SLois Curfman McInnes } 534f1af5d2fSBarry Smith } else if (isdraw) { 535b0a32e0cSBarry Smith PetscDraw draw; 536f1af5d2fSBarry Smith PetscTruth isnull; 537f1af5d2fSBarry Smith 538b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 539b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 540f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 541f1af5d2fSBarry Smith } 54277ed5343SBarry Smith 5438965ea79SLois Curfman McInnes if (size == 1) { 54439ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5453a40ed3dSBarry Smith } else { 5468965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5478965ea79SLois Curfman McInnes Mat A; 548273d9f13SBarry Smith int M = mat->M,N = mat->N,m,row,i,nz,*cols; 54987828ca2SBarry Smith PetscScalar *vals; 5508965ea79SLois Curfman McInnes 5518965ea79SLois Curfman McInnes if (!rank) { 552f1af5d2fSBarry Smith ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 5533a40ed3dSBarry Smith } else { 554f1af5d2fSBarry Smith ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 5558965ea79SLois Curfman McInnes } 556b0a32e0cSBarry Smith PetscLogObjectParent(mat,A); 5578965ea79SLois Curfman McInnes 55839ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 55939ddd567SLois Curfman McInnes but it's quick for now */ 560273d9f13SBarry Smith row = mdn->rstart; m = mdn->A->m; 5618965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 56239ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 56339ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 56439ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 56539ddd567SLois Curfman McInnes row++; 5668965ea79SLois Curfman McInnes } 5678965ea79SLois Curfman McInnes 5686d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5696d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 570b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 571b9b97703SBarry Smith if (!rank) { 5726831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 5738965ea79SLois Curfman McInnes } 574b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 575b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5768965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 5778965ea79SLois Curfman McInnes } 5783a40ed3dSBarry Smith PetscFunctionReturn(0); 5798965ea79SLois Curfman McInnes } 5808965ea79SLois Curfman McInnes 5814a2ae208SSatish Balay #undef __FUNCT__ 5824a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 583b0a32e0cSBarry Smith int MatView_MPIDense(Mat mat,PetscViewer viewer) 5848965ea79SLois Curfman McInnes { 58539ddd567SLois Curfman McInnes int ierr; 586f1af5d2fSBarry Smith PetscTruth isascii,isbinary,isdraw,issocket; 5878965ea79SLois Curfman McInnes 588433994e6SBarry Smith PetscFunctionBegin; 5890f5bd95cSBarry Smith 590b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 591fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 592b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 593fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 5940f5bd95cSBarry Smith 595f1af5d2fSBarry Smith if (isascii || issocket || isdraw) { 596f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 5970f5bd95cSBarry Smith } else if (isbinary) { 5983a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 5995cd90555SBarry Smith } else { 60029bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6018965ea79SLois Curfman McInnes } 6023a40ed3dSBarry Smith PetscFunctionReturn(0); 6038965ea79SLois Curfman McInnes } 6048965ea79SLois Curfman McInnes 6054a2ae208SSatish Balay #undef __FUNCT__ 6064a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 6078f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6088965ea79SLois Curfman McInnes { 6093501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 6103501a2bdSLois Curfman McInnes Mat mdn = mat->A; 6114e220ebcSLois Curfman McInnes int ierr; 612329f5518SBarry Smith PetscReal isend[5],irecv[5]; 6138965ea79SLois Curfman McInnes 6143a40ed3dSBarry Smith PetscFunctionBegin; 615273d9f13SBarry Smith info->rows_global = (double)A->M; 616273d9f13SBarry Smith info->columns_global = (double)A->N; 617273d9f13SBarry Smith info->rows_local = (double)A->m; 618273d9f13SBarry Smith info->columns_local = (double)A->N; 6194e220ebcSLois Curfman McInnes info->block_size = 1.0; 6204e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6214e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6224e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6238965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6244e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6254e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6264e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6274e220ebcSLois Curfman McInnes info->memory = isend[3]; 6284e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6298965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 630d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 6314e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6324e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6334e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6344e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6354e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6368965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 637d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,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 } 6444e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6454e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6464e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6473a40ed3dSBarry Smith PetscFunctionReturn(0); 6488965ea79SLois Curfman McInnes } 6498965ea79SLois Curfman McInnes 6504a2ae208SSatish Balay #undef __FUNCT__ 6514a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 6528f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 6538965ea79SLois Curfman McInnes { 65439ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 655273d9f13SBarry Smith int ierr; 6568965ea79SLois Curfman McInnes 6573a40ed3dSBarry Smith PetscFunctionBegin; 65812c028f9SKris Buschelman switch (op) { 65912c028f9SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 66012c028f9SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 66112c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 66212c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 66312c028f9SKris Buschelman case MAT_COLUMNS_SORTED: 66412c028f9SKris Buschelman case MAT_COLUMNS_UNSORTED: 665273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 66612c028f9SKris Buschelman break; 66712c028f9SKris Buschelman case MAT_ROW_ORIENTED: 668273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 669273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 67012c028f9SKris Buschelman break; 67112c028f9SKris Buschelman case MAT_ROWS_SORTED: 67212c028f9SKris Buschelman case MAT_ROWS_UNSORTED: 67312c028f9SKris Buschelman case MAT_YES_NEW_DIAGONALS: 67412c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 675b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 67612c028f9SKris Buschelman break; 67712c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 678273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 679273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 68012c028f9SKris Buschelman break; 68112c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 682273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 68312c028f9SKris Buschelman break; 68412c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 68529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 686*77e54ba9SKris Buschelman case MAT_SYMMETRIC: 687*77e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 688*77e54ba9SKris Buschelman break; 68912c028f9SKris Buschelman default: 69029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 6913a40ed3dSBarry Smith } 6923a40ed3dSBarry Smith PetscFunctionReturn(0); 6938965ea79SLois Curfman McInnes } 6948965ea79SLois Curfman McInnes 6954a2ae208SSatish Balay #undef __FUNCT__ 6964a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIDense" 69787828ca2SBarry Smith int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v) 6988965ea79SLois Curfman McInnes { 6993501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7003a40ed3dSBarry Smith int lrow,rstart = mat->rstart,rend = mat->rend,ierr; 7018965ea79SLois Curfman McInnes 7023a40ed3dSBarry Smith PetscFunctionBegin; 70329bbc08cSBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 7048965ea79SLois Curfman McInnes lrow = row - rstart; 7053a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 7063a40ed3dSBarry Smith PetscFunctionReturn(0); 7078965ea79SLois Curfman McInnes } 7088965ea79SLois Curfman McInnes 7094a2ae208SSatish Balay #undef __FUNCT__ 7104a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIDense" 71187828ca2SBarry Smith int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v) 7128965ea79SLois Curfman McInnes { 713606d414cSSatish Balay int ierr; 714606d414cSSatish Balay 7153a40ed3dSBarry Smith PetscFunctionBegin; 716606d414cSSatish Balay if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 717606d414cSSatish Balay if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 7183a40ed3dSBarry Smith PetscFunctionReturn(0); 7198965ea79SLois Curfman McInnes } 7208965ea79SLois Curfman McInnes 7214a2ae208SSatish Balay #undef __FUNCT__ 7224a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 7235b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7245b2fa520SLois Curfman McInnes { 7255b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7265b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 72787828ca2SBarry Smith PetscScalar *l,*r,x,*v; 728273d9f13SBarry Smith int ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 7295b2fa520SLois Curfman McInnes 7305b2fa520SLois Curfman McInnes PetscFunctionBegin; 73172d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7325b2fa520SLois Curfman McInnes if (ll) { 73372d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 73429bbc08cSBarry Smith if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 7355b2fa520SLois Curfman McInnes ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7365b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7375b2fa520SLois Curfman McInnes x = l[i]; 7385b2fa520SLois Curfman McInnes v = mat->v + i; 7395b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7405b2fa520SLois Curfman McInnes } 7415b2fa520SLois Curfman McInnes ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 742b0a32e0cSBarry Smith PetscLogFlops(n*m); 7435b2fa520SLois Curfman McInnes } 7445b2fa520SLois Curfman McInnes if (rr) { 74572d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 74629bbc08cSBarry Smith if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 7475b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7485b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7495b2fa520SLois Curfman McInnes ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 7505b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 7515b2fa520SLois Curfman McInnes x = r[i]; 7525b2fa520SLois Curfman McInnes v = mat->v + i*m; 7535b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 7545b2fa520SLois Curfman McInnes } 75572d926a5SLois Curfman McInnes ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 756b0a32e0cSBarry Smith PetscLogFlops(n*m); 7575b2fa520SLois Curfman McInnes } 7585b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 7595b2fa520SLois Curfman McInnes } 7605b2fa520SLois Curfman McInnes 7614a2ae208SSatish Balay #undef __FUNCT__ 7624a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 763064f8208SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 764096963f5SLois Curfman McInnes { 7653501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7663501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 7673501a2bdSLois Curfman McInnes int ierr,i,j; 768329f5518SBarry Smith PetscReal sum = 0.0; 76987828ca2SBarry Smith PetscScalar *v = mat->v; 7703501a2bdSLois Curfman McInnes 7713a40ed3dSBarry Smith PetscFunctionBegin; 7723501a2bdSLois Curfman McInnes if (mdn->size == 1) { 773064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 7743501a2bdSLois Curfman McInnes } else { 7753501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 776273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 777aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 778329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 7793501a2bdSLois Curfman McInnes #else 7803501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 7813501a2bdSLois Curfman McInnes #endif 7823501a2bdSLois Curfman McInnes } 783064f8208SBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 784064f8208SBarry Smith *nrm = sqrt(*nrm); 785b0a32e0cSBarry Smith PetscLogFlops(2*mdn->A->n*mdn->A->m); 7863a40ed3dSBarry Smith } else if (type == NORM_1) { 787329f5518SBarry Smith PetscReal *tmp,*tmp2; 788b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 789273d9f13SBarry Smith tmp2 = tmp + A->N; 790273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 791064f8208SBarry Smith *nrm = 0.0; 7923501a2bdSLois Curfman McInnes v = mat->v; 793273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 794273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 79567e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 7963501a2bdSLois Curfman McInnes } 7973501a2bdSLois Curfman McInnes } 798d7d1e502SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 799273d9f13SBarry Smith for (j=0; j<A->N; j++) { 800064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8013501a2bdSLois Curfman McInnes } 802606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 803b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 8043a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 805329f5518SBarry Smith PetscReal ntemp; 8063501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 807064f8208SBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 8083a40ed3dSBarry Smith } else { 80929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8103501a2bdSLois Curfman McInnes } 8113501a2bdSLois Curfman McInnes } 8123a40ed3dSBarry Smith PetscFunctionReturn(0); 8133501a2bdSLois Curfman McInnes } 8143501a2bdSLois Curfman McInnes 8154a2ae208SSatish Balay #undef __FUNCT__ 8164a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 8178f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 8183501a2bdSLois Curfman McInnes { 8193501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8203501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8213501a2bdSLois Curfman McInnes Mat B; 822273d9f13SBarry Smith int M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 8233501a2bdSLois Curfman McInnes int j,i,ierr; 82487828ca2SBarry Smith PetscScalar *v; 8253501a2bdSLois Curfman McInnes 8263a40ed3dSBarry Smith PetscFunctionBegin; 8277c922b88SBarry Smith if (!matout && M != N) { 82829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 8297056b6fcSBarry Smith } 8307056b6fcSBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 8313501a2bdSLois Curfman McInnes 832273d9f13SBarry Smith m = a->A->m; n = a->A->n; v = Aloc->v; 833b0a32e0cSBarry Smith ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr); 8343501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 8353501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8363501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8373501a2bdSLois Curfman McInnes v += m; 8383501a2bdSLois Curfman McInnes } 839606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 8406d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8416d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8427c922b88SBarry Smith if (matout) { 8433501a2bdSLois Curfman McInnes *matout = B; 8443501a2bdSLois Curfman McInnes } else { 845273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 8463501a2bdSLois Curfman McInnes } 8473a40ed3dSBarry Smith PetscFunctionReturn(0); 848096963f5SLois Curfman McInnes } 849096963f5SLois Curfman McInnes 850d9eff348SSatish Balay #include "petscblaslapack.h" 8514a2ae208SSatish Balay #undef __FUNCT__ 8524a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 853268466fbSBarry Smith int MatScale_MPIDense(const PetscScalar *alpha,Mat inA) 85444cd7ae7SLois Curfman McInnes { 85544cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 85644cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 85744cd7ae7SLois Curfman McInnes int one = 1,nz; 85844cd7ae7SLois Curfman McInnes 8593a40ed3dSBarry Smith PetscFunctionBegin; 860273d9f13SBarry Smith nz = inA->m*inA->N; 861268466fbSBarry Smith BLscal_(&nz,(PetscScalar*)alpha,a->v,&one); 862b0a32e0cSBarry Smith PetscLogFlops(nz); 8633a40ed3dSBarry Smith PetscFunctionReturn(0); 86444cd7ae7SLois Curfman McInnes } 86544cd7ae7SLois Curfman McInnes 8665609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 8678965ea79SLois Curfman McInnes 8684a2ae208SSatish Balay #undef __FUNCT__ 8694a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 870273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A) 871273d9f13SBarry Smith { 872273d9f13SBarry Smith int ierr; 873273d9f13SBarry Smith 874273d9f13SBarry Smith PetscFunctionBegin; 875273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 876273d9f13SBarry Smith PetscFunctionReturn(0); 877273d9f13SBarry Smith } 878273d9f13SBarry Smith 8798965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 88009dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 88109dc0095SBarry Smith MatGetRow_MPIDense, 88209dc0095SBarry Smith MatRestoreRow_MPIDense, 88309dc0095SBarry Smith MatMult_MPIDense, 88497304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 8857c922b88SBarry Smith MatMultTranspose_MPIDense, 8867c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 8878965ea79SLois Curfman McInnes 0, 88809dc0095SBarry Smith 0, 88909dc0095SBarry Smith 0, 89097304618SKris Buschelman /*10*/ 0, 89109dc0095SBarry Smith 0, 89209dc0095SBarry Smith 0, 89309dc0095SBarry Smith 0, 89409dc0095SBarry Smith MatTranspose_MPIDense, 89597304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 89697304618SKris Buschelman 0, 89709dc0095SBarry Smith MatGetDiagonal_MPIDense, 8985b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 89909dc0095SBarry Smith MatNorm_MPIDense, 90097304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 90109dc0095SBarry Smith MatAssemblyEnd_MPIDense, 90209dc0095SBarry Smith 0, 90309dc0095SBarry Smith MatSetOption_MPIDense, 90409dc0095SBarry Smith MatZeroEntries_MPIDense, 90597304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 90609dc0095SBarry Smith 0, 90709dc0095SBarry Smith 0, 90809dc0095SBarry Smith 0, 90909dc0095SBarry Smith 0, 91097304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 911273d9f13SBarry Smith 0, 91209dc0095SBarry Smith 0, 91309dc0095SBarry Smith MatGetArray_MPIDense, 91409dc0095SBarry Smith MatRestoreArray_MPIDense, 91597304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 91609dc0095SBarry Smith 0, 91709dc0095SBarry Smith 0, 91809dc0095SBarry Smith 0, 91909dc0095SBarry Smith 0, 92097304618SKris Buschelman /*40*/ 0, 9212ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 92209dc0095SBarry Smith 0, 92309dc0095SBarry Smith MatGetValues_MPIDense, 92409dc0095SBarry Smith 0, 92597304618SKris Buschelman /*45*/ 0, 92609dc0095SBarry Smith MatScale_MPIDense, 92709dc0095SBarry Smith 0, 92809dc0095SBarry Smith 0, 92909dc0095SBarry Smith 0, 93097304618SKris Buschelman /*50*/ MatGetBlockSize_MPIDense, 93109dc0095SBarry Smith 0, 93209dc0095SBarry Smith 0, 93309dc0095SBarry Smith 0, 93409dc0095SBarry Smith 0, 93597304618SKris Buschelman /*55*/ 0, 93609dc0095SBarry Smith 0, 93709dc0095SBarry Smith 0, 93809dc0095SBarry Smith 0, 93909dc0095SBarry Smith 0, 94097304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 941b9b97703SBarry Smith MatDestroy_MPIDense, 942b9b97703SBarry Smith MatView_MPIDense, 94397304618SKris Buschelman MatGetPetscMaps_Petsc, 94497304618SKris Buschelman 0, 94597304618SKris Buschelman /*65*/ 0, 94697304618SKris Buschelman 0, 94797304618SKris Buschelman 0, 94897304618SKris Buschelman 0, 94997304618SKris Buschelman 0, 95097304618SKris Buschelman /*70*/ 0, 95197304618SKris Buschelman 0, 95297304618SKris Buschelman 0, 95397304618SKris Buschelman 0, 95497304618SKris Buschelman 0, 95597304618SKris Buschelman /*75*/ 0, 95697304618SKris Buschelman 0, 95797304618SKris Buschelman 0, 95897304618SKris Buschelman 0, 95997304618SKris Buschelman 0, 96097304618SKris Buschelman /*80*/ 0, 96197304618SKris Buschelman 0, 96297304618SKris Buschelman 0, 96397304618SKris Buschelman 0, 96497304618SKris Buschelman 0, 96597304618SKris Buschelman /*85*/ MatLoad_MPIDense}; 9668965ea79SLois Curfman McInnes 967273d9f13SBarry Smith EXTERN_C_BEGIN 9684a2ae208SSatish Balay #undef __FUNCT__ 969a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 970a23d5eceSKris Buschelman int MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 971a23d5eceSKris Buschelman { 972a23d5eceSKris Buschelman Mat_MPIDense *a; 973a23d5eceSKris Buschelman int ierr; 974a23d5eceSKris Buschelman 975a23d5eceSKris Buschelman PetscFunctionBegin; 976a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 977a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 978a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 979a23d5eceSKris Buschelman 980a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 981a23d5eceSKris Buschelman ierr = MatCreateSeqDense(PETSC_COMM_SELF,mat->m,mat->N,data,&a->A);CHKERRQ(ierr); 982a23d5eceSKris Buschelman PetscLogObjectParent(mat,a->A); 983a23d5eceSKris Buschelman PetscFunctionReturn(0); 984a23d5eceSKris Buschelman } 985a23d5eceSKris Buschelman EXTERN_C_END 986a23d5eceSKris Buschelman 9870bad9183SKris Buschelman /*MC 988fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 9890bad9183SKris Buschelman 9900bad9183SKris Buschelman Options Database Keys: 9910bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 9920bad9183SKris Buschelman 9930bad9183SKris Buschelman Level: beginner 9940bad9183SKris Buschelman 9950bad9183SKris Buschelman .seealso: MatCreateMPIDense 9960bad9183SKris Buschelman M*/ 9970bad9183SKris Buschelman 998a23d5eceSKris Buschelman EXTERN_C_BEGIN 999a23d5eceSKris Buschelman #undef __FUNCT__ 10004a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1001273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat) 1002273d9f13SBarry Smith { 1003273d9f13SBarry Smith Mat_MPIDense *a; 1004273d9f13SBarry Smith int ierr,i; 1005273d9f13SBarry Smith 1006273d9f13SBarry Smith PetscFunctionBegin; 1007b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1008b0a32e0cSBarry Smith mat->data = (void*)a; 1009273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1010273d9f13SBarry Smith mat->factor = 0; 1011273d9f13SBarry Smith mat->mapping = 0; 1012273d9f13SBarry Smith 1013273d9f13SBarry Smith a->factor = 0; 1014273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1015273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1016273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1017273d9f13SBarry Smith 1018273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1019273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1020273d9f13SBarry Smith a->nvec = mat->n; 1021273d9f13SBarry Smith 1022273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 1023273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 10248a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 10258a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1026273d9f13SBarry Smith 1027273d9f13SBarry Smith /* build local table of row and column ownerships */ 1028b0a32e0cSBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1029273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 1030b0a32e0cSBarry Smith PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1031273d9f13SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1032273d9f13SBarry Smith a->rowners[0] = 0; 1033273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1034273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1035273d9f13SBarry Smith } 1036273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1037273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 1038273d9f13SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1039273d9f13SBarry Smith a->cowners[0] = 0; 1040273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1041273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1042273d9f13SBarry Smith } 1043273d9f13SBarry Smith 1044273d9f13SBarry Smith /* build cache for off array entries formed */ 1045273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1046273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1047273d9f13SBarry Smith 1048273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1049273d9f13SBarry Smith a->lvec = 0; 1050273d9f13SBarry Smith a->Mvctx = 0; 1051273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1052273d9f13SBarry Smith 1053273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1054273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1055273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1056a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1057a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1058a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1059273d9f13SBarry Smith PetscFunctionReturn(0); 1060273d9f13SBarry Smith } 1061273d9f13SBarry Smith EXTERN_C_END 1062273d9f13SBarry Smith 1063209238afSKris Buschelman /*MC 1064002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1065209238afSKris Buschelman 1066209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1067209238afSKris Buschelman and MATMPIDENSE otherwise. 1068209238afSKris Buschelman 1069209238afSKris Buschelman Options Database Keys: 1070209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1071209238afSKris Buschelman 1072209238afSKris Buschelman Level: beginner 1073209238afSKris Buschelman 1074209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1075209238afSKris Buschelman M*/ 1076209238afSKris Buschelman 1077209238afSKris Buschelman EXTERN_C_BEGIN 1078209238afSKris Buschelman #undef __FUNCT__ 1079209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1080209238afSKris Buschelman int MatCreate_Dense(Mat A) { 1081209238afSKris Buschelman int ierr,size; 1082209238afSKris Buschelman 1083209238afSKris Buschelman PetscFunctionBegin; 1084209238afSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1085209238afSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1086209238afSKris Buschelman if (size == 1) { 1087209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1088209238afSKris Buschelman } else { 1089209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1090209238afSKris Buschelman } 1091209238afSKris Buschelman PetscFunctionReturn(0); 1092209238afSKris Buschelman } 1093209238afSKris Buschelman EXTERN_C_END 1094209238afSKris Buschelman 10954a2ae208SSatish Balay #undef __FUNCT__ 10964a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1097273d9f13SBarry Smith /*@C 1098273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1099273d9f13SBarry Smith 1100273d9f13SBarry Smith Not collective 1101273d9f13SBarry Smith 1102273d9f13SBarry Smith Input Parameters: 1103273d9f13SBarry Smith . A - the matrix 1104273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1105273d9f13SBarry Smith to control all matrix memory allocation. 1106273d9f13SBarry Smith 1107273d9f13SBarry Smith Notes: 1108273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1109273d9f13SBarry Smith storage by columns. 1110273d9f13SBarry Smith 1111273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1112273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1113273d9f13SBarry Smith set data=PETSC_NULL. 1114273d9f13SBarry Smith 1115273d9f13SBarry Smith Level: intermediate 1116273d9f13SBarry Smith 1117273d9f13SBarry Smith .keywords: matrix,dense, parallel 1118273d9f13SBarry Smith 1119273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1120273d9f13SBarry Smith @*/ 112187828ca2SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1122273d9f13SBarry Smith { 1123a23d5eceSKris Buschelman int ierr,(*f)(Mat,PetscScalar *); 1124273d9f13SBarry Smith 1125273d9f13SBarry Smith PetscFunctionBegin; 1126565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1127a23d5eceSKris Buschelman if (f) { 1128a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1129a23d5eceSKris Buschelman } 1130273d9f13SBarry Smith PetscFunctionReturn(0); 1131273d9f13SBarry Smith } 1132273d9f13SBarry Smith 11334a2ae208SSatish Balay #undef __FUNCT__ 11344a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 11358965ea79SLois Curfman McInnes /*@C 113639ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 11378965ea79SLois Curfman McInnes 1138db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1139db81eaa0SLois Curfman McInnes 11408965ea79SLois Curfman McInnes Input Parameters: 1141db81eaa0SLois Curfman McInnes + comm - MPI communicator 11428965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1143db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 11448965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1145db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 11467f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1147dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 11488965ea79SLois Curfman McInnes 11498965ea79SLois Curfman McInnes Output Parameter: 1150477f1c0bSLois Curfman McInnes . A - the matrix 11518965ea79SLois Curfman McInnes 1152b259b22eSLois Curfman McInnes Notes: 115339ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 115439ddd567SLois Curfman McInnes storage by columns. 11558965ea79SLois Curfman McInnes 115618f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 115718f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 11587f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 115918f449edSLois Curfman McInnes 11608965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 11618965ea79SLois Curfman McInnes (possibly both). 11628965ea79SLois Curfman McInnes 1163027ccd11SLois Curfman McInnes Level: intermediate 1164027ccd11SLois Curfman McInnes 116539ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 11668965ea79SLois Curfman McInnes 116739ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 11688965ea79SLois Curfman McInnes @*/ 116987828ca2SBarry Smith int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A) 11708965ea79SLois Curfman McInnes { 1171273d9f13SBarry Smith int ierr,size; 11728965ea79SLois Curfman McInnes 11733a40ed3dSBarry Smith PetscFunctionBegin; 1174273d9f13SBarry Smith ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr); 1175273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1176273d9f13SBarry Smith if (size > 1) { 1177273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1178273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1179273d9f13SBarry Smith } else { 1180273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1181273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 11828c469469SLois Curfman McInnes } 11833a40ed3dSBarry Smith PetscFunctionReturn(0); 11848965ea79SLois Curfman McInnes } 11858965ea79SLois Curfman McInnes 11864a2ae208SSatish Balay #undef __FUNCT__ 11874a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 11885609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 11898965ea79SLois Curfman McInnes { 11908965ea79SLois Curfman McInnes Mat mat; 11913501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 119239ddd567SLois Curfman McInnes int ierr; 11938965ea79SLois Curfman McInnes 11943a40ed3dSBarry Smith PetscFunctionBegin; 11958965ea79SLois Curfman McInnes *newmat = 0; 1196273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr); 1197273d9f13SBarry Smith ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr); 1198b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1199b0a32e0cSBarry Smith mat->data = (void*)a; 1200549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 12013501a2bdSLois Curfman McInnes mat->factor = A->factor; 1202c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1203273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 12048965ea79SLois Curfman McInnes 12058965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 12068965ea79SLois Curfman McInnes a->rend = oldmat->rend; 12078965ea79SLois Curfman McInnes a->size = oldmat->size; 12088965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1209e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1210b9b97703SBarry Smith a->nvec = oldmat->nvec; 12113782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1212b0a32e0cSBarry Smith ierr = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1213b0a32e0cSBarry Smith PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1214549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 12158798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 12168965ea79SLois Curfman McInnes 1217329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 12185609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 1219b0a32e0cSBarry Smith PetscLogObjectParent(mat,a->A); 12208965ea79SLois Curfman McInnes *newmat = mat; 12213a40ed3dSBarry Smith PetscFunctionReturn(0); 12228965ea79SLois Curfman McInnes } 12238965ea79SLois Curfman McInnes 1224e090d566SSatish Balay #include "petscsys.h" 12258965ea79SLois Curfman McInnes 12264a2ae208SSatish Balay #undef __FUNCT__ 12274a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 122890ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat) 122990ace30eSBarry Smith { 123040011551SBarry Smith int *rowners,i,size,rank,m,ierr,nz,j; 123187828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 123290ace30eSBarry Smith MPI_Status status; 123390ace30eSBarry Smith 12343a40ed3dSBarry Smith PetscFunctionBegin; 1235d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1236d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 123790ace30eSBarry Smith 123890ace30eSBarry Smith /* determine ownership of all rows */ 123990ace30eSBarry Smith m = M/size + ((M % size) > rank); 1240b0a32e0cSBarry Smith ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1241ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 124290ace30eSBarry Smith rowners[0] = 0; 124390ace30eSBarry Smith for (i=2; i<=size; i++) { 124490ace30eSBarry Smith rowners[i] += rowners[i-1]; 124590ace30eSBarry Smith } 124690ace30eSBarry Smith 124790ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 124890ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 124990ace30eSBarry Smith 125090ace30eSBarry Smith if (!rank) { 125187828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 125290ace30eSBarry Smith 125390ace30eSBarry Smith /* read in my part of the matrix numerical values */ 12540752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 125590ace30eSBarry Smith 125690ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 125790ace30eSBarry Smith vals_ptr = vals; 125890ace30eSBarry Smith for (i=0; i<m; i++) { 125990ace30eSBarry Smith for (j=0; j<N; j++) { 126090ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 126190ace30eSBarry Smith } 126290ace30eSBarry Smith } 126390ace30eSBarry Smith 126490ace30eSBarry Smith /* read in other processors and ship out */ 126590ace30eSBarry Smith for (i=1; i<size; i++) { 126690ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 12670752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1268ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 126990ace30eSBarry Smith } 12703a40ed3dSBarry Smith } else { 127190ace30eSBarry Smith /* receive numeric values */ 127287828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 127390ace30eSBarry Smith 127490ace30eSBarry Smith /* receive message of values*/ 1275ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 127690ace30eSBarry Smith 127790ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 127890ace30eSBarry Smith vals_ptr = vals; 127990ace30eSBarry Smith for (i=0; i<m; i++) { 128090ace30eSBarry Smith for (j=0; j<N; j++) { 128190ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 128290ace30eSBarry Smith } 128390ace30eSBarry Smith } 128490ace30eSBarry Smith } 1285606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1286606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 12876d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12886d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12893a40ed3dSBarry Smith PetscFunctionReturn(0); 129090ace30eSBarry Smith } 129190ace30eSBarry Smith 12924a2ae208SSatish Balay #undef __FUNCT__ 12934a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1294b0a32e0cSBarry Smith int MatLoad_MPIDense(PetscViewer viewer,MatType type,Mat *newmat) 12958965ea79SLois Curfman McInnes { 12968965ea79SLois Curfman McInnes Mat A; 129787828ca2SBarry Smith PetscScalar *vals,*svals; 129819bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 12998965ea79SLois Curfman McInnes MPI_Status status; 13008965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 13018965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols; 130219bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 13033a40ed3dSBarry Smith int i,nz,ierr,j,rstart,rend,fd; 13048965ea79SLois Curfman McInnes 13053a40ed3dSBarry Smith PetscFunctionBegin; 1306d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1307d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13088965ea79SLois Curfman McInnes if (!rank) { 1309b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 13100752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1311552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 13128965ea79SLois Curfman McInnes } 13138965ea79SLois Curfman McInnes 1314ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 131590ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 131690ace30eSBarry Smith 131790ace30eSBarry Smith /* 131890ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 131990ace30eSBarry Smith */ 132090ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 13213a40ed3dSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 13223a40ed3dSBarry Smith PetscFunctionReturn(0); 132390ace30eSBarry Smith } 132490ace30eSBarry Smith 13258965ea79SLois Curfman McInnes /* determine ownership of all rows */ 13268965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 1327b0a32e0cSBarry Smith ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1328ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 13298965ea79SLois Curfman McInnes rowners[0] = 0; 13308965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 13318965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 13328965ea79SLois Curfman McInnes } 13338965ea79SLois Curfman McInnes rstart = rowners[rank]; 13348965ea79SLois Curfman McInnes rend = rowners[rank+1]; 13358965ea79SLois Curfman McInnes 13368965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 1337b0a32e0cSBarry Smith ierr = PetscMalloc(2*(rend-rstart)*sizeof(int),&ourlens);CHKERRQ(ierr); 13388965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 13398965ea79SLois Curfman McInnes if (!rank) { 1340b0a32e0cSBarry Smith ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr); 13410752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 1342b0a32e0cSBarry Smith ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr); 13438965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 1344ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1345606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1346ca161407SBarry Smith } else { 1347ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 13488965ea79SLois Curfman McInnes } 13498965ea79SLois Curfman McInnes 13508965ea79SLois Curfman McInnes if (!rank) { 13518965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 1352b0a32e0cSBarry Smith ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr); 1353549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 13548965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 13558965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 13568965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 13578965ea79SLois Curfman McInnes } 13588965ea79SLois Curfman McInnes } 1359606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 13608965ea79SLois Curfman McInnes 13618965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 13628965ea79SLois Curfman McInnes maxnz = 0; 13638965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 13640452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 13658965ea79SLois Curfman McInnes } 1366b0a32e0cSBarry Smith ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr); 13678965ea79SLois Curfman McInnes 13688965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 13698965ea79SLois Curfman McInnes nz = procsnz[0]; 1370b0a32e0cSBarry Smith ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr); 13710752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 13728965ea79SLois Curfman McInnes 13738965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 13748965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 13758965ea79SLois Curfman McInnes nz = procsnz[i]; 13760752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1377ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 13788965ea79SLois Curfman McInnes } 1379606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 13803a40ed3dSBarry Smith } else { 13818965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 13828965ea79SLois Curfman McInnes nz = 0; 13838965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13848965ea79SLois Curfman McInnes nz += ourlens[i]; 13858965ea79SLois Curfman McInnes } 1386b0a32e0cSBarry Smith ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr); 13878965ea79SLois Curfman McInnes 13888965ea79SLois Curfman McInnes /* receive message of column indices*/ 1389ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1390ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 139129bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 13928965ea79SLois Curfman McInnes } 13938965ea79SLois Curfman McInnes 13948965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1395549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 13968965ea79SLois Curfman McInnes jj = 0; 13978965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 13988965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 13998965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 14008965ea79SLois Curfman McInnes jj++; 14018965ea79SLois Curfman McInnes } 14028965ea79SLois Curfman McInnes } 14038965ea79SLois Curfman McInnes 14048965ea79SLois Curfman McInnes /* create our matrix */ 14058965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14068965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 14078965ea79SLois Curfman McInnes } 1408b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 14098965ea79SLois Curfman McInnes A = *newmat; 14108965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14118965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 14128965ea79SLois Curfman McInnes } 14138965ea79SLois Curfman McInnes 14148965ea79SLois Curfman McInnes if (!rank) { 141587828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14168965ea79SLois Curfman McInnes 14178965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 14188965ea79SLois Curfman McInnes nz = procsnz[0]; 14190752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 14208965ea79SLois Curfman McInnes 14218965ea79SLois Curfman McInnes /* insert into matrix */ 14228965ea79SLois Curfman McInnes jj = rstart; 14238965ea79SLois Curfman McInnes smycols = mycols; 14248965ea79SLois Curfman McInnes svals = vals; 14258965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14268965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14278965ea79SLois Curfman McInnes smycols += ourlens[i]; 14288965ea79SLois Curfman McInnes svals += ourlens[i]; 14298965ea79SLois Curfman McInnes jj++; 14308965ea79SLois Curfman McInnes } 14318965ea79SLois Curfman McInnes 14328965ea79SLois Curfman McInnes /* read in other processors and ship out */ 14338965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14348965ea79SLois Curfman McInnes nz = procsnz[i]; 14350752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1436ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 14378965ea79SLois Curfman McInnes } 1438606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 14393a40ed3dSBarry Smith } else { 14408965ea79SLois Curfman McInnes /* receive numeric values */ 144187828ca2SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14428965ea79SLois Curfman McInnes 14438965ea79SLois Curfman McInnes /* receive message of values*/ 1444ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1445ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 144629bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14478965ea79SLois Curfman McInnes 14488965ea79SLois Curfman McInnes /* insert into matrix */ 14498965ea79SLois Curfman McInnes jj = rstart; 14508965ea79SLois Curfman McInnes smycols = mycols; 14518965ea79SLois Curfman McInnes svals = vals; 14528965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14538965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14548965ea79SLois Curfman McInnes smycols += ourlens[i]; 14558965ea79SLois Curfman McInnes svals += ourlens[i]; 14568965ea79SLois Curfman McInnes jj++; 14578965ea79SLois Curfman McInnes } 14588965ea79SLois Curfman McInnes } 1459606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1460606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1461606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1462606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 14638965ea79SLois Curfman McInnes 14646d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14656d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14663a40ed3dSBarry Smith PetscFunctionReturn(0); 14678965ea79SLois Curfman McInnes } 146890ace30eSBarry Smith 146990ace30eSBarry Smith 147090ace30eSBarry Smith 147190ace30eSBarry Smith 147290ace30eSBarry Smith 1473