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 7*ba8c8a56SBarry Smith #undef __FUNCT__ 8*ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense" 9*ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 10*ba8c8a56SBarry Smith { 11*ba8c8a56SBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 12*ba8c8a56SBarry Smith PetscErrorCode ierr; 13*ba8c8a56SBarry Smith PetscInt lrow,rstart = mat->rstart,rend = mat->rend; 14*ba8c8a56SBarry Smith 15*ba8c8a56SBarry Smith PetscFunctionBegin; 16*ba8c8a56SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 17*ba8c8a56SBarry Smith lrow = row - rstart; 18*ba8c8a56SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr); 19*ba8c8a56SBarry Smith PetscFunctionReturn(0); 20*ba8c8a56SBarry Smith } 21*ba8c8a56SBarry Smith 22*ba8c8a56SBarry Smith #undef __FUNCT__ 23*ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense" 24*ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 25*ba8c8a56SBarry Smith { 26*ba8c8a56SBarry Smith PetscErrorCode ierr; 27*ba8c8a56SBarry Smith 28*ba8c8a56SBarry Smith PetscFunctionBegin; 29*ba8c8a56SBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 30*ba8c8a56SBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 31*ba8c8a56SBarry Smith PetscFunctionReturn(0); 32*ba8c8a56SBarry Smith } 33*ba8c8a56SBarry Smith 340de54da6SSatish Balay EXTERN_C_BEGIN 354a2ae208SSatish Balay #undef __FUNCT__ 364a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 37dfbe8321SBarry Smith PetscErrorCode MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 380de54da6SSatish Balay { 390de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 406849ba73SBarry Smith PetscErrorCode ierr; 4113f74950SBarry Smith PetscInt m = A->m,rstart = mdn->rstart; 4287828ca2SBarry Smith PetscScalar *array; 430de54da6SSatish Balay MPI_Comm comm; 440de54da6SSatish Balay 450de54da6SSatish Balay PetscFunctionBegin; 46273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 470de54da6SSatish Balay 480de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 490de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 500de54da6SSatish Balay 510de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 520de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 535c5985e7SKris Buschelman ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr); 545c5985e7SKris Buschelman ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr); 555c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr); 560de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 570de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 580de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 590de54da6SSatish Balay 600de54da6SSatish Balay *iscopy = PETSC_TRUE; 610de54da6SSatish Balay PetscFunctionReturn(0); 620de54da6SSatish Balay } 630de54da6SSatish Balay EXTERN_C_END 640de54da6SSatish Balay 654a2ae208SSatish Balay #undef __FUNCT__ 664a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 6713f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) 688965ea79SLois Curfman McInnes { 6939b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 70dfbe8321SBarry Smith PetscErrorCode ierr; 7113f74950SBarry Smith PetscInt i,j,rstart = A->rstart,rend = A->rend,row; 72273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 738965ea79SLois Curfman McInnes 743a40ed3dSBarry Smith PetscFunctionBegin; 758965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 765ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 77273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 788965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 798965ea79SLois Curfman McInnes row = idxm[i] - rstart; 8039b7565bSBarry Smith if (roworiented) { 8139b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 823a40ed3dSBarry Smith } else { 838965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 845ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 85273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 8639b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 8739b7565bSBarry Smith } 888965ea79SLois Curfman McInnes } 893a40ed3dSBarry Smith } else { 903782ba37SSatish Balay if (!A->donotstash) { 9139b7565bSBarry Smith if (roworiented) { 928798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 93d36fbae8SSatish Balay } else { 948798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 9539b7565bSBarry Smith } 96b49de8d1SLois Curfman McInnes } 97b49de8d1SLois Curfman McInnes } 983782ba37SSatish Balay } 993a40ed3dSBarry Smith PetscFunctionReturn(0); 100b49de8d1SLois Curfman McInnes } 101b49de8d1SLois Curfman McInnes 1024a2ae208SSatish Balay #undef __FUNCT__ 1034a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 10413f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 105b49de8d1SLois Curfman McInnes { 106b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 107dfbe8321SBarry Smith PetscErrorCode ierr; 10813f74950SBarry Smith PetscInt i,j,rstart = mdn->rstart,rend = mdn->rend,row; 109b49de8d1SLois Curfman McInnes 1103a40ed3dSBarry Smith PetscFunctionBegin; 111b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 11229bbc08cSBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 113273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 114b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 115b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 116b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 11729bbc08cSBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 118273d9f13SBarry Smith if (idxn[j] >= mat->N) { 11929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 120a8c6a408SBarry Smith } 121b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 122b49de8d1SLois Curfman McInnes } 123a8c6a408SBarry Smith } else { 12429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 1258965ea79SLois Curfman McInnes } 1268965ea79SLois Curfman McInnes } 1273a40ed3dSBarry Smith PetscFunctionReturn(0); 1288965ea79SLois Curfman McInnes } 1298965ea79SLois Curfman McInnes 1304a2ae208SSatish Balay #undef __FUNCT__ 1314a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 132dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 133ff14e315SSatish Balay { 134ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 135dfbe8321SBarry Smith PetscErrorCode ierr; 136ff14e315SSatish Balay 1373a40ed3dSBarry Smith PetscFunctionBegin; 138ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1393a40ed3dSBarry Smith PetscFunctionReturn(0); 140ff14e315SSatish Balay } 141ff14e315SSatish Balay 1424a2ae208SSatish Balay #undef __FUNCT__ 1434a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 14413f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 145ca3fa75bSLois Curfman McInnes { 146ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 147ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 1486849ba73SBarry Smith PetscErrorCode ierr; 14913f74950SBarry Smith PetscInt i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 15087828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 151ca3fa75bSLois Curfman McInnes Mat newmat; 152ca3fa75bSLois Curfman McInnes 153ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 154ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 155ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 156b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 157b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 158ca3fa75bSLois Curfman McInnes 159ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1607eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1617eba5e9cSLois Curfman McInnes original matrix! */ 162ca3fa75bSLois Curfman McInnes 163ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1647eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 165ca3fa75bSLois Curfman McInnes 166ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 167ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 16829bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 1697eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 170ca3fa75bSLois Curfman McInnes newmat = *B; 171ca3fa75bSLois Curfman McInnes } else { 172ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 173878740d9SKris Buschelman ierr = MatCreate(A->comm,nrows,cs,PETSC_DECIDE,ncols,&newmat);CHKERRQ(ierr); 174878740d9SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 175878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 176ca3fa75bSLois Curfman McInnes } 177ca3fa75bSLois Curfman McInnes 178ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 179ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 180ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 181ca3fa75bSLois Curfman McInnes 182ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 183ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 184ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 1857eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 186ca3fa75bSLois Curfman McInnes } 187ca3fa75bSLois Curfman McInnes } 188ca3fa75bSLois Curfman McInnes 189ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 190ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 191ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 192ca3fa75bSLois Curfman McInnes 193ca3fa75bSLois Curfman McInnes /* Free work space */ 194ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 195ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 196ca3fa75bSLois Curfman McInnes *B = newmat; 197ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 198ca3fa75bSLois Curfman McInnes } 199ca3fa75bSLois Curfman McInnes 2004a2ae208SSatish Balay #undef __FUNCT__ 2014a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 202dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 203ff14e315SSatish Balay { 2043a40ed3dSBarry Smith PetscFunctionBegin; 2053a40ed3dSBarry Smith PetscFunctionReturn(0); 206ff14e315SSatish Balay } 207ff14e315SSatish Balay 2084a2ae208SSatish Balay #undef __FUNCT__ 2094a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 210dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2118965ea79SLois Curfman McInnes { 21239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 2138965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 214dfbe8321SBarry Smith PetscErrorCode ierr; 21513f74950SBarry Smith PetscInt nstash,reallocs; 2168965ea79SLois Curfman McInnes InsertMode addv; 2178965ea79SLois Curfman McInnes 2183a40ed3dSBarry Smith PetscFunctionBegin; 2198965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 220ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 2217056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 22229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 2238965ea79SLois Curfman McInnes } 224e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2258965ea79SLois Curfman McInnes 2268798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 2278798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 22877431f27SBarry Smith PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs); 2293a40ed3dSBarry Smith PetscFunctionReturn(0); 2308965ea79SLois Curfman McInnes } 2318965ea79SLois Curfman McInnes 2324a2ae208SSatish Balay #undef __FUNCT__ 2334a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 234dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2358965ea79SLois Curfman McInnes { 23639ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2376849ba73SBarry Smith PetscErrorCode ierr; 23813f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 23913f74950SBarry Smith PetscMPIInt n; 24087828ca2SBarry Smith PetscScalar *val; 241e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2428965ea79SLois Curfman McInnes 2433a40ed3dSBarry Smith PetscFunctionBegin; 2448965ea79SLois Curfman McInnes /* wait on receives */ 2457ef1d9bdSSatish Balay while (1) { 2468798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2477ef1d9bdSSatish Balay if (!flg) break; 2488965ea79SLois Curfman McInnes 2497ef1d9bdSSatish Balay for (i=0; i<n;) { 2507ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2517ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2527ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2537ef1d9bdSSatish Balay else ncols = n-i; 2547ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2557ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2567ef1d9bdSSatish Balay i = j; 2578965ea79SLois Curfman McInnes } 2587ef1d9bdSSatish Balay } 2598798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2608965ea79SLois Curfman McInnes 26139ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 26239ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2638965ea79SLois Curfman McInnes 2646d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 26539ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2668965ea79SLois Curfman McInnes } 2673a40ed3dSBarry Smith PetscFunctionReturn(0); 2688965ea79SLois Curfman McInnes } 2698965ea79SLois Curfman McInnes 2704a2ae208SSatish Balay #undef __FUNCT__ 2714a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 272dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 2738965ea79SLois Curfman McInnes { 274dfbe8321SBarry Smith PetscErrorCode ierr; 27539ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2763a40ed3dSBarry Smith 2773a40ed3dSBarry Smith PetscFunctionBegin; 2783a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2793a40ed3dSBarry Smith PetscFunctionReturn(0); 2808965ea79SLois Curfman McInnes } 2818965ea79SLois Curfman McInnes 2828965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2838965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2848965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2853501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2868965ea79SLois Curfman McInnes routine. 2878965ea79SLois Curfman McInnes */ 2884a2ae208SSatish Balay #undef __FUNCT__ 2894a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 290dfbe8321SBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag) 2918965ea79SLois Curfman McInnes { 29239ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2936849ba73SBarry Smith PetscErrorCode ierr; 29413f74950SBarry Smith PetscInt i,N,*rows,*owners = l->rowners; 29513f74950SBarry Smith PetscInt *nprocs,j,idx,nsends; 29613f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 29713f74950SBarry Smith PetscInt *rvalues,tag = A->tag,count,base,slen,*source; 29813f74950SBarry Smith PetscInt *lens,*lrows,*values; 29913f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 3008965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 3018965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3028965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 3038965ea79SLois Curfman McInnes IS istmp; 30435d8aa7fSBarry Smith PetscTruth found; 3058965ea79SLois Curfman McInnes 3063a40ed3dSBarry Smith PetscFunctionBegin; 307b9b97703SBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 3088965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 3098965ea79SLois Curfman McInnes 3108965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 31113f74950SBarry Smith ierr = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr); 31213f74950SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 31313f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/ 3148965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3158965ea79SLois Curfman McInnes idx = rows[i]; 31635d8aa7fSBarry Smith found = PETSC_FALSE; 3178965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3188965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 319c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3208965ea79SLois Curfman McInnes } 3218965ea79SLois Curfman McInnes } 32229bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3238965ea79SLois Curfman McInnes } 324c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 3258965ea79SLois Curfman McInnes 3268965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 327c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3288965ea79SLois Curfman McInnes 3298965ea79SLois Curfman McInnes /* post receives: */ 33013f74950SBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr); 331b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3328965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 33313f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3348965ea79SLois Curfman McInnes } 3358965ea79SLois Curfman McInnes 3368965ea79SLois Curfman McInnes /* do sends: 3378965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3388965ea79SLois Curfman McInnes the ith processor 3398965ea79SLois Curfman McInnes */ 34013f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr); 341b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 34213f74950SBarry Smith ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr); 3438965ea79SLois Curfman McInnes starts[0] = 0; 344c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3458965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3468965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3478965ea79SLois Curfman McInnes } 3488965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3498965ea79SLois Curfman McInnes 3508965ea79SLois Curfman McInnes starts[0] = 0; 351c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3528965ea79SLois Curfman McInnes count = 0; 3538965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 354c1dc657dSBarry Smith if (nprocs[2*i+1]) { 35513f74950SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3568965ea79SLois Curfman McInnes } 3578965ea79SLois Curfman McInnes } 358606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3598965ea79SLois Curfman McInnes 3608965ea79SLois Curfman McInnes base = owners[rank]; 3618965ea79SLois Curfman McInnes 3628965ea79SLois Curfman McInnes /* wait on receives */ 36313f74950SBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr); 3648965ea79SLois Curfman McInnes source = lens + nrecvs; 3658965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3668965ea79SLois Curfman McInnes while (count) { 367ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3688965ea79SLois Curfman McInnes /* unpack receives into our local space */ 36913f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 3708965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3718965ea79SLois Curfman McInnes lens[imdex] = n; 3728965ea79SLois Curfman McInnes slen += n; 3738965ea79SLois Curfman McInnes count--; 3748965ea79SLois Curfman McInnes } 375606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3768965ea79SLois Curfman McInnes 3778965ea79SLois Curfman McInnes /* move the data into the send scatter */ 37813f74950SBarry Smith ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr); 3798965ea79SLois Curfman McInnes count = 0; 3808965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 3818965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3828965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 3838965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3848965ea79SLois Curfman McInnes } 3858965ea79SLois Curfman McInnes } 386606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 387606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 388606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 389606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3908965ea79SLois Curfman McInnes 3918965ea79SLois Curfman McInnes /* actually zap the local rows */ 392029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 393b0a32e0cSBarry Smith PetscLogObjectParent(A,istmp); 394606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3958965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 3968965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 3978965ea79SLois Curfman McInnes 3988965ea79SLois Curfman McInnes /* wait on sends */ 3998965ea79SLois Curfman McInnes if (nsends) { 400b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 401ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 402606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 4038965ea79SLois Curfman McInnes } 404606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 405606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 4068965ea79SLois Curfman McInnes 4073a40ed3dSBarry Smith PetscFunctionReturn(0); 4088965ea79SLois Curfman McInnes } 4098965ea79SLois Curfman McInnes 4104a2ae208SSatish Balay #undef __FUNCT__ 4114a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 412dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4138965ea79SLois Curfman McInnes { 41439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 415dfbe8321SBarry Smith PetscErrorCode ierr; 416c456f294SBarry Smith 4173a40ed3dSBarry Smith PetscFunctionBegin; 41843a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41943a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 42044cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4213a40ed3dSBarry Smith PetscFunctionReturn(0); 4228965ea79SLois Curfman McInnes } 4238965ea79SLois Curfman McInnes 4244a2ae208SSatish Balay #undef __FUNCT__ 4254a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 426dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4278965ea79SLois Curfman McInnes { 42839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 429dfbe8321SBarry Smith PetscErrorCode ierr; 430c456f294SBarry Smith 4313a40ed3dSBarry Smith PetscFunctionBegin; 43243a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 43343a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 43444cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4353a40ed3dSBarry Smith PetscFunctionReturn(0); 4368965ea79SLois Curfman McInnes } 4378965ea79SLois Curfman McInnes 4384a2ae208SSatish Balay #undef __FUNCT__ 4394a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 440dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 441096963f5SLois Curfman McInnes { 442096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 443dfbe8321SBarry Smith PetscErrorCode ierr; 44487828ca2SBarry Smith PetscScalar zero = 0.0; 445096963f5SLois Curfman McInnes 4463a40ed3dSBarry Smith PetscFunctionBegin; 4473501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 4487c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 449537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 450537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4513a40ed3dSBarry Smith PetscFunctionReturn(0); 452096963f5SLois Curfman McInnes } 453096963f5SLois Curfman McInnes 4544a2ae208SSatish Balay #undef __FUNCT__ 4554a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 456dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 457096963f5SLois Curfman McInnes { 458096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 459dfbe8321SBarry Smith PetscErrorCode ierr; 460096963f5SLois Curfman McInnes 4613a40ed3dSBarry Smith PetscFunctionBegin; 4623501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4637c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 464537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 465537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4663a40ed3dSBarry Smith PetscFunctionReturn(0); 467096963f5SLois Curfman McInnes } 468096963f5SLois Curfman McInnes 4694a2ae208SSatish Balay #undef __FUNCT__ 4704a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 471dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 4728965ea79SLois Curfman McInnes { 47339ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 474096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 475dfbe8321SBarry Smith PetscErrorCode ierr; 47613f74950SBarry Smith PetscInt len,i,n,m = A->m,radd; 47787828ca2SBarry Smith PetscScalar *x,zero = 0.0; 478ed3cc1f0SBarry Smith 4793a40ed3dSBarry Smith PetscFunctionBegin; 480273d9f13SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 4811ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 482096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 483273d9f13SBarry Smith if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 484273d9f13SBarry Smith len = PetscMin(a->A->m,a->A->n); 4857ddc982cSLois Curfman McInnes radd = a->rstart*m; 48644cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 487096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 488096963f5SLois Curfman McInnes } 4891ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 4903a40ed3dSBarry Smith PetscFunctionReturn(0); 4918965ea79SLois Curfman McInnes } 4928965ea79SLois Curfman McInnes 4934a2ae208SSatish Balay #undef __FUNCT__ 4944a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 495dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 4968965ea79SLois Curfman McInnes { 4973501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 498dfbe8321SBarry Smith PetscErrorCode ierr; 499ed3cc1f0SBarry Smith 5003a40ed3dSBarry Smith PetscFunctionBegin; 50194d884c6SBarry Smith 502aa482453SBarry Smith #if defined(PETSC_USE_LOG) 50377431f27SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N); 5048965ea79SLois Curfman McInnes #endif 5058798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 506606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 5073501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 5083501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 5093501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 510622d7880SLois Curfman McInnes if (mdn->factor) { 511606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 512606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 513606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 514606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 515622d7880SLois Curfman McInnes } 516606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 517901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr); 518901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 5193a40ed3dSBarry Smith PetscFunctionReturn(0); 5208965ea79SLois Curfman McInnes } 52139ddd567SLois Curfman McInnes 5224a2ae208SSatish Balay #undef __FUNCT__ 5234a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 5246849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5258965ea79SLois Curfman McInnes { 52639ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 527dfbe8321SBarry Smith PetscErrorCode ierr; 5287056b6fcSBarry Smith 5293a40ed3dSBarry Smith PetscFunctionBegin; 53039ddd567SLois Curfman McInnes if (mdn->size == 1) { 53139ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5328965ea79SLois Curfman McInnes } 53329bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5343a40ed3dSBarry Smith PetscFunctionReturn(0); 5358965ea79SLois Curfman McInnes } 5368965ea79SLois Curfman McInnes 5374a2ae208SSatish Balay #undef __FUNCT__ 5384a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 5396849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5408965ea79SLois Curfman McInnes { 54139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 542dfbe8321SBarry Smith PetscErrorCode ierr; 54332dcc486SBarry Smith PetscMPIInt size = mdn->size,rank = mdn->rank; 544b0a32e0cSBarry Smith PetscViewerType vtype; 54532077d6dSBarry Smith PetscTruth iascii,isdraw; 546b0a32e0cSBarry Smith PetscViewer sviewer; 547f3ef73ceSBarry Smith PetscViewerFormat format; 5488965ea79SLois Curfman McInnes 5493a40ed3dSBarry Smith PetscFunctionBegin; 55032077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 551fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 55232077d6dSBarry Smith if (iascii) { 553b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 554b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 555456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 5564e220ebcSLois Curfman McInnes MatInfo info; 557888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 55877431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m, 55977431f27SBarry Smith (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 560b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5613501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5623a40ed3dSBarry Smith PetscFunctionReturn(0); 563fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 5643a40ed3dSBarry Smith PetscFunctionReturn(0); 5658965ea79SLois Curfman McInnes } 566f1af5d2fSBarry Smith } else if (isdraw) { 567b0a32e0cSBarry Smith PetscDraw draw; 568f1af5d2fSBarry Smith PetscTruth isnull; 569f1af5d2fSBarry Smith 570b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 571b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 572f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 573f1af5d2fSBarry Smith } 57477ed5343SBarry Smith 5758965ea79SLois Curfman McInnes if (size == 1) { 57639ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5773a40ed3dSBarry Smith } else { 5788965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5798965ea79SLois Curfman McInnes Mat A; 58013f74950SBarry Smith PetscInt M = mat->M,N = mat->N,m,row,i,nz; 581*ba8c8a56SBarry Smith PetscInt *cols; 582*ba8c8a56SBarry Smith PetscScalar *vals; 5838965ea79SLois Curfman McInnes 5848965ea79SLois Curfman McInnes if (!rank) { 585878740d9SKris Buschelman ierr = MatCreate(mat->comm,M,N,M,N,&A);CHKERRQ(ierr); 5863a40ed3dSBarry Smith } else { 587878740d9SKris Buschelman ierr = MatCreate(mat->comm,0,0,M,N,&A);CHKERRQ(ierr); 5888965ea79SLois Curfman McInnes } 589be5d1d56SKris Buschelman /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */ 590878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 591878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 592b0a32e0cSBarry Smith PetscLogObjectParent(mat,A); 5938965ea79SLois Curfman McInnes 59439ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 59539ddd567SLois Curfman McInnes but it's quick for now */ 596273d9f13SBarry Smith row = mdn->rstart; m = mdn->A->m; 5978965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 598*ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 599*ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 600*ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 60139ddd567SLois Curfman McInnes row++; 6028965ea79SLois Curfman McInnes } 6038965ea79SLois Curfman McInnes 6046d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6056d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 606b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 607b9b97703SBarry Smith if (!rank) { 6086831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 6098965ea79SLois Curfman McInnes } 610b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 611b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6128965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 6138965ea79SLois Curfman McInnes } 6143a40ed3dSBarry Smith PetscFunctionReturn(0); 6158965ea79SLois Curfman McInnes } 6168965ea79SLois Curfman McInnes 6174a2ae208SSatish Balay #undef __FUNCT__ 6184a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 619dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 6208965ea79SLois Curfman McInnes { 621dfbe8321SBarry Smith PetscErrorCode ierr; 62232077d6dSBarry Smith PetscTruth iascii,isbinary,isdraw,issocket; 6238965ea79SLois Curfman McInnes 624433994e6SBarry Smith PetscFunctionBegin; 6250f5bd95cSBarry Smith 62632077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 627fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 628b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 629fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 6300f5bd95cSBarry Smith 63132077d6dSBarry Smith if (iascii || issocket || isdraw) { 632f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 6330f5bd95cSBarry Smith } else if (isbinary) { 6343a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6355cd90555SBarry Smith } else { 636958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6378965ea79SLois Curfman McInnes } 6383a40ed3dSBarry Smith PetscFunctionReturn(0); 6398965ea79SLois Curfman McInnes } 6408965ea79SLois Curfman McInnes 6414a2ae208SSatish Balay #undef __FUNCT__ 6424a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 643dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6448965ea79SLois Curfman McInnes { 6453501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 6463501a2bdSLois Curfman McInnes Mat mdn = mat->A; 647dfbe8321SBarry Smith PetscErrorCode ierr; 648329f5518SBarry Smith PetscReal isend[5],irecv[5]; 6498965ea79SLois Curfman McInnes 6503a40ed3dSBarry Smith PetscFunctionBegin; 651273d9f13SBarry Smith info->rows_global = (double)A->M; 652273d9f13SBarry Smith info->columns_global = (double)A->N; 653273d9f13SBarry Smith info->rows_local = (double)A->m; 654273d9f13SBarry Smith info->columns_local = (double)A->N; 6554e220ebcSLois Curfman McInnes info->block_size = 1.0; 6564e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6574e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6584e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6598965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6604e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6614e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6624e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6634e220ebcSLois Curfman McInnes info->memory = isend[3]; 6644e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6658965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 666d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 6674e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6684e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6694e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6704e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6714e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6728965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 673d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 6744e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6754e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6764e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6774e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6784e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6798965ea79SLois Curfman McInnes } 6804e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6814e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6824e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6833a40ed3dSBarry Smith PetscFunctionReturn(0); 6848965ea79SLois Curfman McInnes } 6858965ea79SLois Curfman McInnes 6864a2ae208SSatish Balay #undef __FUNCT__ 6874a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 688dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op) 6898965ea79SLois Curfman McInnes { 69039ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 691dfbe8321SBarry Smith PetscErrorCode ierr; 6928965ea79SLois Curfman McInnes 6933a40ed3dSBarry Smith PetscFunctionBegin; 69412c028f9SKris Buschelman switch (op) { 69512c028f9SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 69612c028f9SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 69712c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 69812c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 69912c028f9SKris Buschelman case MAT_COLUMNS_SORTED: 70012c028f9SKris Buschelman case MAT_COLUMNS_UNSORTED: 701273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 70212c028f9SKris Buschelman break; 70312c028f9SKris Buschelman case MAT_ROW_ORIENTED: 704273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 705273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 70612c028f9SKris Buschelman break; 70712c028f9SKris Buschelman case MAT_ROWS_SORTED: 70812c028f9SKris Buschelman case MAT_ROWS_UNSORTED: 70912c028f9SKris Buschelman case MAT_YES_NEW_DIAGONALS: 71012c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 711b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 71212c028f9SKris Buschelman break; 71312c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 714273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 715273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 71612c028f9SKris Buschelman break; 71712c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 718273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 71912c028f9SKris Buschelman break; 72012c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 72129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 72277e54ba9SKris Buschelman case MAT_SYMMETRIC: 72377e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 7249a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 7259a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 7269a4540c5SBarry Smith case MAT_HERMITIAN: 7279a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 7289a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 7299a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 73077e54ba9SKris Buschelman break; 73112c028f9SKris Buschelman default: 73229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 7333a40ed3dSBarry Smith } 7343a40ed3dSBarry Smith PetscFunctionReturn(0); 7358965ea79SLois Curfman McInnes } 7368965ea79SLois Curfman McInnes 7378965ea79SLois Curfman McInnes 7384a2ae208SSatish Balay #undef __FUNCT__ 7394a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 740dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7415b2fa520SLois Curfman McInnes { 7425b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7435b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 74487828ca2SBarry Smith PetscScalar *l,*r,x,*v; 745dfbe8321SBarry Smith PetscErrorCode ierr; 74613f74950SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 7475b2fa520SLois Curfman McInnes 7485b2fa520SLois Curfman McInnes PetscFunctionBegin; 74972d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7505b2fa520SLois Curfman McInnes if (ll) { 75172d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 75229bbc08cSBarry Smith if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 7531ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7545b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7555b2fa520SLois Curfman McInnes x = l[i]; 7565b2fa520SLois Curfman McInnes v = mat->v + i; 7575b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7585b2fa520SLois Curfman McInnes } 7591ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 760b0a32e0cSBarry Smith PetscLogFlops(n*m); 7615b2fa520SLois Curfman McInnes } 7625b2fa520SLois Curfman McInnes if (rr) { 76372d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 76429bbc08cSBarry Smith if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 7655b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7665b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7671ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 7685b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 7695b2fa520SLois Curfman McInnes x = r[i]; 7705b2fa520SLois Curfman McInnes v = mat->v + i*m; 7715b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 7725b2fa520SLois Curfman McInnes } 7731ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 774b0a32e0cSBarry Smith PetscLogFlops(n*m); 7755b2fa520SLois Curfman McInnes } 7765b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 7775b2fa520SLois Curfman McInnes } 7785b2fa520SLois Curfman McInnes 7794a2ae208SSatish Balay #undef __FUNCT__ 7804a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 781dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 782096963f5SLois Curfman McInnes { 7833501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7843501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 785dfbe8321SBarry Smith PetscErrorCode ierr; 78613f74950SBarry Smith PetscInt i,j; 787329f5518SBarry Smith PetscReal sum = 0.0; 78887828ca2SBarry Smith PetscScalar *v = mat->v; 7893501a2bdSLois Curfman McInnes 7903a40ed3dSBarry Smith PetscFunctionBegin; 7913501a2bdSLois Curfman McInnes if (mdn->size == 1) { 792064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 7933501a2bdSLois Curfman McInnes } else { 7943501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 795273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 796aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 797329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 7983501a2bdSLois Curfman McInnes #else 7993501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 8003501a2bdSLois Curfman McInnes #endif 8013501a2bdSLois Curfman McInnes } 802064f8208SBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 803064f8208SBarry Smith *nrm = sqrt(*nrm); 804b0a32e0cSBarry Smith PetscLogFlops(2*mdn->A->n*mdn->A->m); 8053a40ed3dSBarry Smith } else if (type == NORM_1) { 806329f5518SBarry Smith PetscReal *tmp,*tmp2; 807b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 808273d9f13SBarry Smith tmp2 = tmp + A->N; 809273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 810064f8208SBarry Smith *nrm = 0.0; 8113501a2bdSLois Curfman McInnes v = mat->v; 812273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 813273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 81467e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8153501a2bdSLois Curfman McInnes } 8163501a2bdSLois Curfman McInnes } 817d7d1e502SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 818273d9f13SBarry Smith for (j=0; j<A->N; j++) { 819064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8203501a2bdSLois Curfman McInnes } 821606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 822b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 8233a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 824329f5518SBarry Smith PetscReal ntemp; 8253501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 826064f8208SBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 8273a40ed3dSBarry Smith } else { 82829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8293501a2bdSLois Curfman McInnes } 8303501a2bdSLois Curfman McInnes } 8313a40ed3dSBarry Smith PetscFunctionReturn(0); 8323501a2bdSLois Curfman McInnes } 8333501a2bdSLois Curfman McInnes 8344a2ae208SSatish Balay #undef __FUNCT__ 8354a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 836dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout) 8373501a2bdSLois Curfman McInnes { 8383501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8393501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8403501a2bdSLois Curfman McInnes Mat B; 84113f74950SBarry Smith PetscInt M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 8426849ba73SBarry Smith PetscErrorCode ierr; 84313f74950SBarry Smith PetscInt j,i; 84487828ca2SBarry Smith PetscScalar *v; 8453501a2bdSLois Curfman McInnes 8463a40ed3dSBarry Smith PetscFunctionBegin; 8477c922b88SBarry Smith if (!matout && M != N) { 84829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 8497056b6fcSBarry Smith } 850878740d9SKris Buschelman ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B);CHKERRQ(ierr); 851878740d9SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 852878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr); 8533501a2bdSLois Curfman McInnes 854273d9f13SBarry Smith m = a->A->m; n = a->A->n; v = Aloc->v; 85513f74950SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr); 8563501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 8573501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8583501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8593501a2bdSLois Curfman McInnes v += m; 8603501a2bdSLois Curfman McInnes } 861606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 8626d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8636d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8647c922b88SBarry Smith if (matout) { 8653501a2bdSLois Curfman McInnes *matout = B; 8663501a2bdSLois Curfman McInnes } else { 867273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 8683501a2bdSLois Curfman McInnes } 8693a40ed3dSBarry Smith PetscFunctionReturn(0); 870096963f5SLois Curfman McInnes } 871096963f5SLois Curfman McInnes 872d9eff348SSatish Balay #include "petscblaslapack.h" 8734a2ae208SSatish Balay #undef __FUNCT__ 8744a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 875dfbe8321SBarry Smith PetscErrorCode MatScale_MPIDense(const PetscScalar *alpha,Mat inA) 87644cd7ae7SLois Curfman McInnes { 87744cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 87844cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 8794ce68768SBarry Smith PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N; 88044cd7ae7SLois Curfman McInnes 8813a40ed3dSBarry Smith PetscFunctionBegin; 882268466fbSBarry Smith BLscal_(&nz,(PetscScalar*)alpha,a->v,&one); 883b0a32e0cSBarry Smith PetscLogFlops(nz); 8843a40ed3dSBarry Smith PetscFunctionReturn(0); 88544cd7ae7SLois Curfman McInnes } 88644cd7ae7SLois Curfman McInnes 8876849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 8888965ea79SLois Curfman McInnes 8894a2ae208SSatish Balay #undef __FUNCT__ 8904a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 891dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A) 892273d9f13SBarry Smith { 893dfbe8321SBarry Smith PetscErrorCode ierr; 894273d9f13SBarry Smith 895273d9f13SBarry Smith PetscFunctionBegin; 896273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 897273d9f13SBarry Smith PetscFunctionReturn(0); 898273d9f13SBarry Smith } 899273d9f13SBarry Smith 9008965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 90109dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 90209dc0095SBarry Smith MatGetRow_MPIDense, 90309dc0095SBarry Smith MatRestoreRow_MPIDense, 90409dc0095SBarry Smith MatMult_MPIDense, 90597304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 9067c922b88SBarry Smith MatMultTranspose_MPIDense, 9077c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 9088965ea79SLois Curfman McInnes 0, 90909dc0095SBarry Smith 0, 91009dc0095SBarry Smith 0, 91197304618SKris Buschelman /*10*/ 0, 91209dc0095SBarry Smith 0, 91309dc0095SBarry Smith 0, 91409dc0095SBarry Smith 0, 91509dc0095SBarry Smith MatTranspose_MPIDense, 91697304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 91797304618SKris Buschelman 0, 91809dc0095SBarry Smith MatGetDiagonal_MPIDense, 9195b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 92009dc0095SBarry Smith MatNorm_MPIDense, 92197304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 92209dc0095SBarry Smith MatAssemblyEnd_MPIDense, 92309dc0095SBarry Smith 0, 92409dc0095SBarry Smith MatSetOption_MPIDense, 92509dc0095SBarry Smith MatZeroEntries_MPIDense, 92697304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 92709dc0095SBarry Smith 0, 92809dc0095SBarry Smith 0, 92909dc0095SBarry Smith 0, 93009dc0095SBarry Smith 0, 93197304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 932273d9f13SBarry Smith 0, 93309dc0095SBarry Smith 0, 93409dc0095SBarry Smith MatGetArray_MPIDense, 93509dc0095SBarry Smith MatRestoreArray_MPIDense, 93697304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 93709dc0095SBarry Smith 0, 93809dc0095SBarry Smith 0, 93909dc0095SBarry Smith 0, 94009dc0095SBarry Smith 0, 94197304618SKris Buschelman /*40*/ 0, 9422ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 94309dc0095SBarry Smith 0, 94409dc0095SBarry Smith MatGetValues_MPIDense, 94509dc0095SBarry Smith 0, 94697304618SKris Buschelman /*45*/ 0, 94709dc0095SBarry Smith MatScale_MPIDense, 94809dc0095SBarry Smith 0, 94909dc0095SBarry Smith 0, 95009dc0095SBarry Smith 0, 951521d7252SBarry Smith /*50*/ 0, 95209dc0095SBarry Smith 0, 95309dc0095SBarry Smith 0, 95409dc0095SBarry Smith 0, 95509dc0095SBarry Smith 0, 95697304618SKris Buschelman /*55*/ 0, 95709dc0095SBarry Smith 0, 95809dc0095SBarry Smith 0, 95909dc0095SBarry Smith 0, 96009dc0095SBarry Smith 0, 96197304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 962b9b97703SBarry Smith MatDestroy_MPIDense, 963b9b97703SBarry Smith MatView_MPIDense, 96497304618SKris Buschelman MatGetPetscMaps_Petsc, 96597304618SKris Buschelman 0, 96697304618SKris Buschelman /*65*/ 0, 96797304618SKris Buschelman 0, 96897304618SKris Buschelman 0, 96997304618SKris Buschelman 0, 97097304618SKris Buschelman 0, 97197304618SKris Buschelman /*70*/ 0, 97297304618SKris Buschelman 0, 97397304618SKris Buschelman 0, 97497304618SKris Buschelman 0, 97597304618SKris Buschelman 0, 97697304618SKris Buschelman /*75*/ 0, 97797304618SKris Buschelman 0, 97897304618SKris Buschelman 0, 97997304618SKris Buschelman 0, 98097304618SKris Buschelman 0, 98197304618SKris Buschelman /*80*/ 0, 98297304618SKris Buschelman 0, 98397304618SKris Buschelman 0, 98497304618SKris Buschelman 0, 985865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense, 986865e5f61SKris Buschelman 0, 987865e5f61SKris Buschelman 0, 988865e5f61SKris Buschelman 0, 989865e5f61SKris Buschelman 0, 990865e5f61SKris Buschelman 0, 991865e5f61SKris Buschelman /*90*/ 0, 992865e5f61SKris Buschelman 0, 993865e5f61SKris Buschelman 0, 994865e5f61SKris Buschelman 0, 995865e5f61SKris Buschelman 0, 996865e5f61SKris Buschelman /*95*/ 0, 997865e5f61SKris Buschelman 0, 998865e5f61SKris Buschelman 0, 999865e5f61SKris Buschelman 0}; 10008965ea79SLois Curfman McInnes 1001273d9f13SBarry Smith EXTERN_C_BEGIN 10024a2ae208SSatish Balay #undef __FUNCT__ 1003a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 1004dfbe8321SBarry Smith PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1005a23d5eceSKris Buschelman { 1006a23d5eceSKris Buschelman Mat_MPIDense *a; 1007dfbe8321SBarry Smith PetscErrorCode ierr; 1008a23d5eceSKris Buschelman 1009a23d5eceSKris Buschelman PetscFunctionBegin; 1010a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1011a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1012a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1013a23d5eceSKris Buschelman 1014a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 10155c5985e7SKris Buschelman ierr = MatCreate(PETSC_COMM_SELF,mat->m,mat->N,mat->m,mat->N,&a->A);CHKERRQ(ierr); 10165c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 10175c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 1018a23d5eceSKris Buschelman PetscLogObjectParent(mat,a->A); 1019a23d5eceSKris Buschelman PetscFunctionReturn(0); 1020a23d5eceSKris Buschelman } 1021a23d5eceSKris Buschelman EXTERN_C_END 1022a23d5eceSKris Buschelman 10230bad9183SKris Buschelman /*MC 1024fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 10250bad9183SKris Buschelman 10260bad9183SKris Buschelman Options Database Keys: 10270bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 10280bad9183SKris Buschelman 10290bad9183SKris Buschelman Level: beginner 10300bad9183SKris Buschelman 10310bad9183SKris Buschelman .seealso: MatCreateMPIDense 10320bad9183SKris Buschelman M*/ 10330bad9183SKris Buschelman 1034a23d5eceSKris Buschelman EXTERN_C_BEGIN 1035a23d5eceSKris Buschelman #undef __FUNCT__ 10364a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1037dfbe8321SBarry Smith PetscErrorCode MatCreate_MPIDense(Mat mat) 1038273d9f13SBarry Smith { 1039273d9f13SBarry Smith Mat_MPIDense *a; 1040dfbe8321SBarry Smith PetscErrorCode ierr; 104113f74950SBarry Smith PetscInt i; 1042273d9f13SBarry Smith 1043273d9f13SBarry Smith PetscFunctionBegin; 1044b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1045b0a32e0cSBarry Smith mat->data = (void*)a; 1046273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1047273d9f13SBarry Smith mat->factor = 0; 1048273d9f13SBarry Smith mat->mapping = 0; 1049273d9f13SBarry Smith 1050273d9f13SBarry Smith a->factor = 0; 1051273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1052273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1053273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1054273d9f13SBarry Smith 1055273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1056273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1057273d9f13SBarry Smith a->nvec = mat->n; 1058273d9f13SBarry Smith 1059273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 1060273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 10618a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 10628a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1063273d9f13SBarry Smith 1064273d9f13SBarry Smith /* build local table of row and column ownerships */ 106513f74950SBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 1066273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 106713f74950SBarry Smith PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 106813f74950SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1069273d9f13SBarry Smith a->rowners[0] = 0; 1070273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1071273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1072273d9f13SBarry Smith } 1073273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1074273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 107513f74950SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1076273d9f13SBarry Smith a->cowners[0] = 0; 1077273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1078273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1079273d9f13SBarry Smith } 1080273d9f13SBarry Smith 1081273d9f13SBarry Smith /* build cache for off array entries formed */ 1082273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1083273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1084273d9f13SBarry Smith 1085273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1086273d9f13SBarry Smith a->lvec = 0; 1087273d9f13SBarry Smith a->Mvctx = 0; 1088273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1089273d9f13SBarry Smith 1090273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1091273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1092273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1093a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1094a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1095a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1096273d9f13SBarry Smith PetscFunctionReturn(0); 1097273d9f13SBarry Smith } 1098273d9f13SBarry Smith EXTERN_C_END 1099273d9f13SBarry Smith 1100209238afSKris Buschelman /*MC 1101002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1102209238afSKris Buschelman 1103209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1104209238afSKris Buschelman and MATMPIDENSE otherwise. 1105209238afSKris Buschelman 1106209238afSKris Buschelman Options Database Keys: 1107209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1108209238afSKris Buschelman 1109209238afSKris Buschelman Level: beginner 1110209238afSKris Buschelman 1111209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1112209238afSKris Buschelman M*/ 1113209238afSKris Buschelman 1114209238afSKris Buschelman EXTERN_C_BEGIN 1115209238afSKris Buschelman #undef __FUNCT__ 1116209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1117dfbe8321SBarry Smith PetscErrorCode MatCreate_Dense(Mat A) 1118dfbe8321SBarry Smith { 11196849ba73SBarry Smith PetscErrorCode ierr; 112013f74950SBarry Smith PetscMPIInt size; 1121209238afSKris Buschelman 1122209238afSKris Buschelman PetscFunctionBegin; 1123209238afSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1124209238afSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1125209238afSKris Buschelman if (size == 1) { 1126209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1127209238afSKris Buschelman } else { 1128209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1129209238afSKris Buschelman } 1130209238afSKris Buschelman PetscFunctionReturn(0); 1131209238afSKris Buschelman } 1132209238afSKris Buschelman EXTERN_C_END 1133209238afSKris Buschelman 11344a2ae208SSatish Balay #undef __FUNCT__ 11354a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1136273d9f13SBarry Smith /*@C 1137273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1138273d9f13SBarry Smith 1139273d9f13SBarry Smith Not collective 1140273d9f13SBarry Smith 1141273d9f13SBarry Smith Input Parameters: 1142273d9f13SBarry Smith . A - the matrix 1143273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1144273d9f13SBarry Smith to control all matrix memory allocation. 1145273d9f13SBarry Smith 1146273d9f13SBarry Smith Notes: 1147273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1148273d9f13SBarry Smith storage by columns. 1149273d9f13SBarry Smith 1150273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1151273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1152273d9f13SBarry Smith set data=PETSC_NULL. 1153273d9f13SBarry Smith 1154273d9f13SBarry Smith Level: intermediate 1155273d9f13SBarry Smith 1156273d9f13SBarry Smith .keywords: matrix,dense, parallel 1157273d9f13SBarry Smith 1158273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1159273d9f13SBarry Smith @*/ 1160dfbe8321SBarry Smith PetscErrorCode MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1161273d9f13SBarry Smith { 1162dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar *); 1163273d9f13SBarry Smith 1164273d9f13SBarry Smith PetscFunctionBegin; 1165565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1166a23d5eceSKris Buschelman if (f) { 1167a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1168a23d5eceSKris Buschelman } 1169273d9f13SBarry Smith PetscFunctionReturn(0); 1170273d9f13SBarry Smith } 1171273d9f13SBarry Smith 11724a2ae208SSatish Balay #undef __FUNCT__ 11734a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 11748965ea79SLois Curfman McInnes /*@C 117539ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 11768965ea79SLois Curfman McInnes 1177db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1178db81eaa0SLois Curfman McInnes 11798965ea79SLois Curfman McInnes Input Parameters: 1180db81eaa0SLois Curfman McInnes + comm - MPI communicator 11818965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1182db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 11838965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1184db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 11857f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1186dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 11878965ea79SLois Curfman McInnes 11888965ea79SLois Curfman McInnes Output Parameter: 1189477f1c0bSLois Curfman McInnes . A - the matrix 11908965ea79SLois Curfman McInnes 1191b259b22eSLois Curfman McInnes Notes: 119239ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 119339ddd567SLois Curfman McInnes storage by columns. 11948965ea79SLois Curfman McInnes 119518f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 119618f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 11977f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 119818f449edSLois Curfman McInnes 11998965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 12008965ea79SLois Curfman McInnes (possibly both). 12018965ea79SLois Curfman McInnes 1202027ccd11SLois Curfman McInnes Level: intermediate 1203027ccd11SLois Curfman McInnes 120439ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 12058965ea79SLois Curfman McInnes 120639ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 12078965ea79SLois Curfman McInnes @*/ 120813f74950SBarry Smith PetscErrorCode MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 12098965ea79SLois Curfman McInnes { 12106849ba73SBarry Smith PetscErrorCode ierr; 121113f74950SBarry Smith PetscMPIInt size; 12128965ea79SLois Curfman McInnes 12133a40ed3dSBarry Smith PetscFunctionBegin; 1214273d9f13SBarry Smith ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr); 1215273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1216273d9f13SBarry Smith if (size > 1) { 1217273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1218273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1219273d9f13SBarry Smith } else { 1220273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1221273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 12228c469469SLois Curfman McInnes } 12233a40ed3dSBarry Smith PetscFunctionReturn(0); 12248965ea79SLois Curfman McInnes } 12258965ea79SLois Curfman McInnes 12264a2ae208SSatish Balay #undef __FUNCT__ 12274a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 12286849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 12298965ea79SLois Curfman McInnes { 12308965ea79SLois Curfman McInnes Mat mat; 12313501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1232dfbe8321SBarry Smith PetscErrorCode ierr; 12338965ea79SLois Curfman McInnes 12343a40ed3dSBarry Smith PetscFunctionBegin; 12358965ea79SLois Curfman McInnes *newmat = 0; 1236273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr); 1237be5d1d56SKris Buschelman ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr); 1238b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1239b0a32e0cSBarry Smith mat->data = (void*)a; 1240549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 12413501a2bdSLois Curfman McInnes mat->factor = A->factor; 1242c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1243273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 12448965ea79SLois Curfman McInnes 12458965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 12468965ea79SLois Curfman McInnes a->rend = oldmat->rend; 12478965ea79SLois Curfman McInnes a->size = oldmat->size; 12488965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1249e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1250b9b97703SBarry Smith a->nvec = oldmat->nvec; 12513782ba37SSatish Balay a->donotstash = oldmat->donotstash; 125213f74950SBarry Smith ierr = PetscMalloc((a->size+1)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 125313f74950SBarry Smith PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 125413f74950SBarry Smith ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr); 12558798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 12568965ea79SLois Curfman McInnes 1257329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 12585609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 1259b0a32e0cSBarry Smith PetscLogObjectParent(mat,a->A); 12608965ea79SLois Curfman McInnes *newmat = mat; 12613a40ed3dSBarry Smith PetscFunctionReturn(0); 12628965ea79SLois Curfman McInnes } 12638965ea79SLois Curfman McInnes 1264e090d566SSatish Balay #include "petscsys.h" 12658965ea79SLois Curfman McInnes 12664a2ae208SSatish Balay #undef __FUNCT__ 12674a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 126813f74950SBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,const MatType type,Mat *newmat) 126990ace30eSBarry Smith { 12706849ba73SBarry Smith PetscErrorCode ierr; 127132dcc486SBarry Smith PetscMPIInt rank,size; 127213f74950SBarry Smith PetscInt *rowners,i,m,nz,j; 127387828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 127490ace30eSBarry Smith MPI_Status status; 127590ace30eSBarry Smith 12763a40ed3dSBarry Smith PetscFunctionBegin; 1277d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1278d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 127990ace30eSBarry Smith 128090ace30eSBarry Smith /* determine ownership of all rows */ 128190ace30eSBarry Smith m = M/size + ((M % size) > rank); 128213f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 128313f74950SBarry Smith ierr = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 128490ace30eSBarry Smith rowners[0] = 0; 128590ace30eSBarry Smith for (i=2; i<=size; i++) { 128690ace30eSBarry Smith rowners[i] += rowners[i-1]; 128790ace30eSBarry Smith } 128890ace30eSBarry Smith 1289878740d9SKris Buschelman ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 1290be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1291878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 129290ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 129390ace30eSBarry Smith 129490ace30eSBarry Smith if (!rank) { 129587828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 129690ace30eSBarry Smith 129790ace30eSBarry Smith /* read in my part of the matrix numerical values */ 12980752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 129990ace30eSBarry Smith 130090ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 130190ace30eSBarry Smith vals_ptr = vals; 130290ace30eSBarry Smith for (i=0; i<m; i++) { 130390ace30eSBarry Smith for (j=0; j<N; j++) { 130490ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 130590ace30eSBarry Smith } 130690ace30eSBarry Smith } 130790ace30eSBarry Smith 130890ace30eSBarry Smith /* read in other processors and ship out */ 130990ace30eSBarry Smith for (i=1; i<size; i++) { 131090ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 13110752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1312ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 131390ace30eSBarry Smith } 13143a40ed3dSBarry Smith } else { 131590ace30eSBarry Smith /* receive numeric values */ 131687828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 131790ace30eSBarry Smith 131890ace30eSBarry Smith /* receive message of values*/ 1319ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 132090ace30eSBarry Smith 132190ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 132290ace30eSBarry Smith vals_ptr = vals; 132390ace30eSBarry Smith for (i=0; i<m; i++) { 132490ace30eSBarry Smith for (j=0; j<N; j++) { 132590ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 132690ace30eSBarry Smith } 132790ace30eSBarry Smith } 132890ace30eSBarry Smith } 1329606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1330606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 13316d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13326d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13333a40ed3dSBarry Smith PetscFunctionReturn(0); 133490ace30eSBarry Smith } 133590ace30eSBarry Smith 13364a2ae208SSatish Balay #undef __FUNCT__ 13374a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1338dfbe8321SBarry Smith PetscErrorCode MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat) 13398965ea79SLois Curfman McInnes { 13408965ea79SLois Curfman McInnes Mat A; 134187828ca2SBarry Smith PetscScalar *vals,*svals; 134219bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 13438965ea79SLois Curfman McInnes MPI_Status status; 134413f74950SBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz; 134513f74950SBarry Smith PetscInt header[4],*rowlengths = 0,M,N,*cols; 134613f74950SBarry Smith PetscInt *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols; 134713f74950SBarry Smith PetscInt i,nz,j,rstart,rend; 134813f74950SBarry Smith int fd; 13496849ba73SBarry Smith PetscErrorCode ierr; 13508965ea79SLois Curfman McInnes 13513a40ed3dSBarry Smith PetscFunctionBegin; 1352d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1353d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13548965ea79SLois Curfman McInnes if (!rank) { 1355b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 13560752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1357552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 13588965ea79SLois Curfman McInnes } 13598965ea79SLois Curfman McInnes 136013f74950SBarry Smith ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 136190ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 136290ace30eSBarry Smith 136390ace30eSBarry Smith /* 136490ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 136590ace30eSBarry Smith */ 136690ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1367be5d1d56SKris Buschelman ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 13683a40ed3dSBarry Smith PetscFunctionReturn(0); 136990ace30eSBarry Smith } 137090ace30eSBarry Smith 13718965ea79SLois Curfman McInnes /* determine ownership of all rows */ 13728965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 137313f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 1374ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 13758965ea79SLois Curfman McInnes rowners[0] = 0; 13768965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 13778965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 13788965ea79SLois Curfman McInnes } 13798965ea79SLois Curfman McInnes rstart = rowners[rank]; 13808965ea79SLois Curfman McInnes rend = rowners[rank+1]; 13818965ea79SLois Curfman McInnes 13828965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 138313f74950SBarry Smith ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr); 13848965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 13858965ea79SLois Curfman McInnes if (!rank) { 138613f74950SBarry Smith ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 13870752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 138813f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr); 13898965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 13901466f1e1SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 1391606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1392ca161407SBarry Smith } else { 13931466f1e1SBarry Smith ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 13948965ea79SLois Curfman McInnes } 13958965ea79SLois Curfman McInnes 13968965ea79SLois Curfman McInnes if (!rank) { 13978965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 139813f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr); 139913f74950SBarry Smith ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 14008965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14018965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 14028965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 14038965ea79SLois Curfman McInnes } 14048965ea79SLois Curfman McInnes } 1405606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 14068965ea79SLois Curfman McInnes 14078965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 14088965ea79SLois Curfman McInnes maxnz = 0; 14098965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14100452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 14118965ea79SLois Curfman McInnes } 141213f74950SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr); 14138965ea79SLois Curfman McInnes 14148965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 14158965ea79SLois Curfman McInnes nz = procsnz[0]; 141613f74950SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14170752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 14188965ea79SLois Curfman McInnes 14198965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 14208965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14218965ea79SLois Curfman McInnes nz = procsnz[i]; 14220752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 142313f74950SBarry Smith ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 14248965ea79SLois Curfman McInnes } 1425606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 14263a40ed3dSBarry Smith } else { 14278965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 14288965ea79SLois Curfman McInnes nz = 0; 14298965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14308965ea79SLois Curfman McInnes nz += ourlens[i]; 14318965ea79SLois Curfman McInnes } 143213f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14338965ea79SLois Curfman McInnes 14348965ea79SLois Curfman McInnes /* receive message of column indices*/ 143513f74950SBarry Smith ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 143613f74950SBarry Smith ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 143729bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14388965ea79SLois Curfman McInnes } 14398965ea79SLois Curfman McInnes 14408965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 144113f74950SBarry Smith ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr); 14428965ea79SLois Curfman McInnes jj = 0; 14438965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14448965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 14458965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 14468965ea79SLois Curfman McInnes jj++; 14478965ea79SLois Curfman McInnes } 14488965ea79SLois Curfman McInnes } 14498965ea79SLois Curfman McInnes 14508965ea79SLois Curfman McInnes /* create our matrix */ 14518965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14528965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 14538965ea79SLois Curfman McInnes } 1454878740d9SKris Buschelman ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 1455be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1456878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 14578965ea79SLois Curfman McInnes A = *newmat; 14588965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14598965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 14608965ea79SLois Curfman McInnes } 14618965ea79SLois Curfman McInnes 14628965ea79SLois Curfman McInnes if (!rank) { 146387828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14648965ea79SLois Curfman McInnes 14658965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 14668965ea79SLois Curfman McInnes nz = procsnz[0]; 14670752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 14688965ea79SLois Curfman McInnes 14698965ea79SLois Curfman McInnes /* insert into matrix */ 14708965ea79SLois Curfman McInnes jj = rstart; 14718965ea79SLois Curfman McInnes smycols = mycols; 14728965ea79SLois Curfman McInnes svals = vals; 14738965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14748965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14758965ea79SLois Curfman McInnes smycols += ourlens[i]; 14768965ea79SLois Curfman McInnes svals += ourlens[i]; 14778965ea79SLois Curfman McInnes jj++; 14788965ea79SLois Curfman McInnes } 14798965ea79SLois Curfman McInnes 14808965ea79SLois Curfman McInnes /* read in other processors and ship out */ 14818965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14828965ea79SLois Curfman McInnes nz = procsnz[i]; 14830752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1484ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 14858965ea79SLois Curfman McInnes } 1486606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 14873a40ed3dSBarry Smith } else { 14888965ea79SLois Curfman McInnes /* receive numeric values */ 148984d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14908965ea79SLois Curfman McInnes 14918965ea79SLois Curfman McInnes /* receive message of values*/ 1492ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1493ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 149429bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14958965ea79SLois Curfman McInnes 14968965ea79SLois Curfman McInnes /* insert into matrix */ 14978965ea79SLois Curfman McInnes jj = rstart; 14988965ea79SLois Curfman McInnes smycols = mycols; 14998965ea79SLois Curfman McInnes svals = vals; 15008965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15018965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 15028965ea79SLois Curfman McInnes smycols += ourlens[i]; 15038965ea79SLois Curfman McInnes svals += ourlens[i]; 15048965ea79SLois Curfman McInnes jj++; 15058965ea79SLois Curfman McInnes } 15068965ea79SLois Curfman McInnes } 1507606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1508606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1509606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1510606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 15118965ea79SLois Curfman McInnes 15126d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15136d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15143a40ed3dSBarry Smith PetscFunctionReturn(0); 15158965ea79SLois Curfman McInnes } 151690ace30eSBarry Smith 151790ace30eSBarry Smith 151890ace30eSBarry Smith 151990ace30eSBarry Smith 152090ace30eSBarry Smith 1521