1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 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" 88965ea79SLois Curfman McInnes 9ba8c8a56SBarry Smith #undef __FUNCT__ 10ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense" 11ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 12ba8c8a56SBarry Smith { 13ba8c8a56SBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 14ba8c8a56SBarry Smith PetscErrorCode ierr; 15ba8c8a56SBarry Smith PetscInt lrow,rstart = mat->rstart,rend = mat->rend; 16ba8c8a56SBarry Smith 17ba8c8a56SBarry Smith PetscFunctionBegin; 18ba8c8a56SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 19ba8c8a56SBarry Smith lrow = row - rstart; 20ba8c8a56SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr); 21ba8c8a56SBarry Smith PetscFunctionReturn(0); 22ba8c8a56SBarry Smith } 23ba8c8a56SBarry Smith 24ba8c8a56SBarry Smith #undef __FUNCT__ 25ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense" 26ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 27ba8c8a56SBarry Smith { 28ba8c8a56SBarry Smith PetscErrorCode ierr; 29ba8c8a56SBarry Smith 30ba8c8a56SBarry Smith PetscFunctionBegin; 31ba8c8a56SBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 32ba8c8a56SBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 33ba8c8a56SBarry Smith PetscFunctionReturn(0); 34ba8c8a56SBarry Smith } 35ba8c8a56SBarry Smith 360de54da6SSatish Balay EXTERN_C_BEGIN 374a2ae208SSatish Balay #undef __FUNCT__ 384a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 39be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 400de54da6SSatish Balay { 410de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 426849ba73SBarry Smith PetscErrorCode ierr; 4313f74950SBarry Smith PetscInt m = A->m,rstart = mdn->rstart; 4487828ca2SBarry Smith PetscScalar *array; 450de54da6SSatish Balay MPI_Comm comm; 460de54da6SSatish Balay 470de54da6SSatish Balay PetscFunctionBegin; 48273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 490de54da6SSatish Balay 500de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 510de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 520de54da6SSatish Balay 530de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 540de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 55*f69a0ea3SMatthew Knepley ierr = MatCreate(comm,B);CHKERRQ(ierr); 56*f69a0ea3SMatthew Knepley ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr); 575c5985e7SKris Buschelman ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr); 585c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr); 590de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 600de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 610de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 620de54da6SSatish Balay 630de54da6SSatish Balay *iscopy = PETSC_TRUE; 640de54da6SSatish Balay PetscFunctionReturn(0); 650de54da6SSatish Balay } 660de54da6SSatish Balay EXTERN_C_END 670de54da6SSatish Balay 684a2ae208SSatish Balay #undef __FUNCT__ 694a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 7013f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) 718965ea79SLois Curfman McInnes { 7239b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 73dfbe8321SBarry Smith PetscErrorCode ierr; 7413f74950SBarry Smith PetscInt i,j,rstart = A->rstart,rend = A->rend,row; 75273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 768965ea79SLois Curfman McInnes 773a40ed3dSBarry Smith PetscFunctionBegin; 788965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 795ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 80273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 818965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 828965ea79SLois Curfman McInnes row = idxm[i] - rstart; 8339b7565bSBarry Smith if (roworiented) { 8439b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 853a40ed3dSBarry Smith } else { 868965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 875ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 88273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 8939b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 9039b7565bSBarry Smith } 918965ea79SLois Curfman McInnes } 923a40ed3dSBarry Smith } else { 933782ba37SSatish Balay if (!A->donotstash) { 9439b7565bSBarry Smith if (roworiented) { 958798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 96d36fbae8SSatish Balay } else { 978798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 9839b7565bSBarry Smith } 99b49de8d1SLois Curfman McInnes } 100b49de8d1SLois Curfman McInnes } 1013782ba37SSatish Balay } 1023a40ed3dSBarry Smith PetscFunctionReturn(0); 103b49de8d1SLois Curfman McInnes } 104b49de8d1SLois Curfman McInnes 1054a2ae208SSatish Balay #undef __FUNCT__ 1064a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 10713f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 108b49de8d1SLois Curfman McInnes { 109b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 110dfbe8321SBarry Smith PetscErrorCode ierr; 11113f74950SBarry Smith PetscInt i,j,rstart = mdn->rstart,rend = mdn->rend,row; 112b49de8d1SLois Curfman McInnes 1133a40ed3dSBarry Smith PetscFunctionBegin; 114b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 11529bbc08cSBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 116273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 117b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 118b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 119b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 12029bbc08cSBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 121273d9f13SBarry Smith if (idxn[j] >= mat->N) { 12229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 123a8c6a408SBarry Smith } 124b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 125b49de8d1SLois Curfman McInnes } 126a8c6a408SBarry Smith } else { 12729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 1288965ea79SLois Curfman McInnes } 1298965ea79SLois Curfman McInnes } 1303a40ed3dSBarry Smith PetscFunctionReturn(0); 1318965ea79SLois Curfman McInnes } 1328965ea79SLois Curfman McInnes 1334a2ae208SSatish Balay #undef __FUNCT__ 1344a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 135dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 136ff14e315SSatish Balay { 137ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 138dfbe8321SBarry Smith PetscErrorCode ierr; 139ff14e315SSatish Balay 1403a40ed3dSBarry Smith PetscFunctionBegin; 141ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1423a40ed3dSBarry Smith PetscFunctionReturn(0); 143ff14e315SSatish Balay } 144ff14e315SSatish Balay 1454a2ae208SSatish Balay #undef __FUNCT__ 1464a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 14713f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 148ca3fa75bSLois Curfman McInnes { 149ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 150ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 1516849ba73SBarry Smith PetscErrorCode ierr; 15213f74950SBarry Smith PetscInt i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 15387828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 154ca3fa75bSLois Curfman McInnes Mat newmat; 155ca3fa75bSLois Curfman McInnes 156ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 157ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 158ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 159b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 160b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 161ca3fa75bSLois Curfman McInnes 162ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1637eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1647eba5e9cSLois Curfman McInnes original matrix! */ 165ca3fa75bSLois Curfman McInnes 166ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1677eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 168ca3fa75bSLois Curfman McInnes 169ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 170ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 17129bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 1727eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 173ca3fa75bSLois Curfman McInnes newmat = *B; 174ca3fa75bSLois Curfman McInnes } else { 175ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 176*f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr); 177*f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr); 178878740d9SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 179878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 180ca3fa75bSLois Curfman McInnes } 181ca3fa75bSLois Curfman McInnes 182ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 183ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 184ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 185ca3fa75bSLois Curfman McInnes 186ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 187ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 188ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 1897eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 190ca3fa75bSLois Curfman McInnes } 191ca3fa75bSLois Curfman McInnes } 192ca3fa75bSLois Curfman McInnes 193ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 194ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 195ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 196ca3fa75bSLois Curfman McInnes 197ca3fa75bSLois Curfman McInnes /* Free work space */ 198ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 199ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 200ca3fa75bSLois Curfman McInnes *B = newmat; 201ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 202ca3fa75bSLois Curfman McInnes } 203ca3fa75bSLois Curfman McInnes 2044a2ae208SSatish Balay #undef __FUNCT__ 2054a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 206dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 207ff14e315SSatish Balay { 2083a40ed3dSBarry Smith PetscFunctionBegin; 2093a40ed3dSBarry Smith PetscFunctionReturn(0); 210ff14e315SSatish Balay } 211ff14e315SSatish Balay 2124a2ae208SSatish Balay #undef __FUNCT__ 2134a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 214dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2158965ea79SLois Curfman McInnes { 21639ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 2178965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 218dfbe8321SBarry Smith PetscErrorCode ierr; 21913f74950SBarry Smith PetscInt nstash,reallocs; 2208965ea79SLois Curfman McInnes InsertMode addv; 2218965ea79SLois Curfman McInnes 2223a40ed3dSBarry Smith PetscFunctionBegin; 2238965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 224ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 2257056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 22629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 2278965ea79SLois Curfman McInnes } 228e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2298965ea79SLois Curfman McInnes 2308798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 2318798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 23263ba0a88SBarry Smith ierr = PetscLogInfo((mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs));CHKERRQ(ierr); 2333a40ed3dSBarry Smith PetscFunctionReturn(0); 2348965ea79SLois Curfman McInnes } 2358965ea79SLois Curfman McInnes 2364a2ae208SSatish Balay #undef __FUNCT__ 2374a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 238dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2398965ea79SLois Curfman McInnes { 24039ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2416849ba73SBarry Smith PetscErrorCode ierr; 24213f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 24313f74950SBarry Smith PetscMPIInt n; 24487828ca2SBarry Smith PetscScalar *val; 245e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2468965ea79SLois Curfman McInnes 2473a40ed3dSBarry Smith PetscFunctionBegin; 2488965ea79SLois Curfman McInnes /* wait on receives */ 2497ef1d9bdSSatish Balay while (1) { 2508798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2517ef1d9bdSSatish Balay if (!flg) break; 2528965ea79SLois Curfman McInnes 2537ef1d9bdSSatish Balay for (i=0; i<n;) { 2547ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2557ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2567ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2577ef1d9bdSSatish Balay else ncols = n-i; 2587ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2597ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2607ef1d9bdSSatish Balay i = j; 2618965ea79SLois Curfman McInnes } 2627ef1d9bdSSatish Balay } 2638798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2648965ea79SLois Curfman McInnes 26539ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 26639ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2678965ea79SLois Curfman McInnes 2686d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 26939ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2708965ea79SLois Curfman McInnes } 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 2728965ea79SLois Curfman McInnes } 2738965ea79SLois Curfman McInnes 2744a2ae208SSatish Balay #undef __FUNCT__ 2754a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 276dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 2778965ea79SLois Curfman McInnes { 278dfbe8321SBarry Smith PetscErrorCode ierr; 27939ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2803a40ed3dSBarry Smith 2813a40ed3dSBarry Smith PetscFunctionBegin; 2823a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2833a40ed3dSBarry Smith PetscFunctionReturn(0); 2848965ea79SLois Curfman McInnes } 2858965ea79SLois Curfman McInnes 2868965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2878965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2888965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2893501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2908965ea79SLois Curfman McInnes routine. 2918965ea79SLois Curfman McInnes */ 2924a2ae208SSatish Balay #undef __FUNCT__ 2934a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 294dfbe8321SBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag) 2958965ea79SLois Curfman McInnes { 29639ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2976849ba73SBarry Smith PetscErrorCode ierr; 29813f74950SBarry Smith PetscInt i,N,*rows,*owners = l->rowners; 29913f74950SBarry Smith PetscInt *nprocs,j,idx,nsends; 30013f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 30113f74950SBarry Smith PetscInt *rvalues,tag = A->tag,count,base,slen,*source; 30213f74950SBarry Smith PetscInt *lens,*lrows,*values; 30313f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 3048965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 3058965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3068965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 3078965ea79SLois Curfman McInnes IS istmp; 30835d8aa7fSBarry Smith PetscTruth found; 3098965ea79SLois Curfman McInnes 3103a40ed3dSBarry Smith PetscFunctionBegin; 311b9b97703SBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 3128965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 3138965ea79SLois Curfman McInnes 3148965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 31513f74950SBarry Smith ierr = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr); 31613f74950SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 31713f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/ 3188965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3198965ea79SLois Curfman McInnes idx = rows[i]; 32035d8aa7fSBarry Smith found = PETSC_FALSE; 3218965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3228965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 323c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3248965ea79SLois Curfman McInnes } 3258965ea79SLois Curfman McInnes } 32629bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3278965ea79SLois Curfman McInnes } 328c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 3298965ea79SLois Curfman McInnes 3308965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 331c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3328965ea79SLois Curfman McInnes 3338965ea79SLois Curfman McInnes /* post receives: */ 33413f74950SBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr); 335b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3368965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 33713f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3388965ea79SLois Curfman McInnes } 3398965ea79SLois Curfman McInnes 3408965ea79SLois Curfman McInnes /* do sends: 3418965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3428965ea79SLois Curfman McInnes the ith processor 3438965ea79SLois Curfman McInnes */ 34413f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr); 345b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 34613f74950SBarry Smith ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr); 3478965ea79SLois Curfman McInnes starts[0] = 0; 348c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3498965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3508965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3518965ea79SLois Curfman McInnes } 3528965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3538965ea79SLois Curfman McInnes 3548965ea79SLois Curfman McInnes starts[0] = 0; 355c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3568965ea79SLois Curfman McInnes count = 0; 3578965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 358c1dc657dSBarry Smith if (nprocs[2*i+1]) { 35913f74950SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3608965ea79SLois Curfman McInnes } 3618965ea79SLois Curfman McInnes } 362606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3638965ea79SLois Curfman McInnes 3648965ea79SLois Curfman McInnes base = owners[rank]; 3658965ea79SLois Curfman McInnes 3668965ea79SLois Curfman McInnes /* wait on receives */ 36713f74950SBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr); 3688965ea79SLois Curfman McInnes source = lens + nrecvs; 3698965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3708965ea79SLois Curfman McInnes while (count) { 371ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3728965ea79SLois Curfman McInnes /* unpack receives into our local space */ 37313f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 3748965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3758965ea79SLois Curfman McInnes lens[imdex] = n; 3768965ea79SLois Curfman McInnes slen += n; 3778965ea79SLois Curfman McInnes count--; 3788965ea79SLois Curfman McInnes } 379606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3808965ea79SLois Curfman McInnes 3818965ea79SLois Curfman McInnes /* move the data into the send scatter */ 38213f74950SBarry Smith ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr); 3838965ea79SLois Curfman McInnes count = 0; 3848965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 3858965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3868965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 3878965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3888965ea79SLois Curfman McInnes } 3898965ea79SLois Curfman McInnes } 390606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 391606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 392606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 393606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3948965ea79SLois Curfman McInnes 3958965ea79SLois Curfman McInnes /* actually zap the local rows */ 396029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 39752e6d16bSBarry Smith ierr = PetscLogObjectParent(A,istmp);CHKERRQ(ierr); 398606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3998965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 4008965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 4018965ea79SLois Curfman McInnes 4028965ea79SLois Curfman McInnes /* wait on sends */ 4038965ea79SLois Curfman McInnes if (nsends) { 404b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 405ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 406606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 4078965ea79SLois Curfman McInnes } 408606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 409606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 4108965ea79SLois Curfman McInnes 4113a40ed3dSBarry Smith PetscFunctionReturn(0); 4128965ea79SLois Curfman McInnes } 4138965ea79SLois Curfman McInnes 4144a2ae208SSatish Balay #undef __FUNCT__ 4154a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 416dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4178965ea79SLois Curfman McInnes { 41839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 419dfbe8321SBarry Smith PetscErrorCode ierr; 420c456f294SBarry Smith 4213a40ed3dSBarry Smith PetscFunctionBegin; 42243a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 42343a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 42444cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4253a40ed3dSBarry Smith PetscFunctionReturn(0); 4268965ea79SLois Curfman McInnes } 4278965ea79SLois Curfman McInnes 4284a2ae208SSatish Balay #undef __FUNCT__ 4294a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 430dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4318965ea79SLois Curfman McInnes { 43239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 433dfbe8321SBarry Smith PetscErrorCode ierr; 434c456f294SBarry Smith 4353a40ed3dSBarry Smith PetscFunctionBegin; 43643a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 43743a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 43844cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4393a40ed3dSBarry Smith PetscFunctionReturn(0); 4408965ea79SLois Curfman McInnes } 4418965ea79SLois Curfman McInnes 4424a2ae208SSatish Balay #undef __FUNCT__ 4434a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 444dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 445096963f5SLois Curfman McInnes { 446096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 447dfbe8321SBarry Smith PetscErrorCode ierr; 44887828ca2SBarry Smith PetscScalar zero = 0.0; 449096963f5SLois Curfman McInnes 4503a40ed3dSBarry Smith PetscFunctionBegin; 4512dcb1b2aSMatthew Knepley ierr = VecSet(yy,zero);CHKERRQ(ierr); 4527c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 453537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 454537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4553a40ed3dSBarry Smith PetscFunctionReturn(0); 456096963f5SLois Curfman McInnes } 457096963f5SLois Curfman McInnes 4584a2ae208SSatish Balay #undef __FUNCT__ 4594a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 460dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 461096963f5SLois Curfman McInnes { 462096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 463dfbe8321SBarry Smith PetscErrorCode ierr; 464096963f5SLois Curfman McInnes 4653a40ed3dSBarry Smith PetscFunctionBegin; 4663501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4677c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 468537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 469537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4703a40ed3dSBarry Smith PetscFunctionReturn(0); 471096963f5SLois Curfman McInnes } 472096963f5SLois Curfman McInnes 4734a2ae208SSatish Balay #undef __FUNCT__ 4744a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 475dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 4768965ea79SLois Curfman McInnes { 47739ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 478096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 479dfbe8321SBarry Smith PetscErrorCode ierr; 48013f74950SBarry Smith PetscInt len,i,n,m = A->m,radd; 48187828ca2SBarry Smith PetscScalar *x,zero = 0.0; 482ed3cc1f0SBarry Smith 4833a40ed3dSBarry Smith PetscFunctionBegin; 4842dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 4851ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 486096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 487273d9f13SBarry Smith if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 488273d9f13SBarry Smith len = PetscMin(a->A->m,a->A->n); 4897ddc982cSLois Curfman McInnes radd = a->rstart*m; 49044cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 491096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 492096963f5SLois Curfman McInnes } 4931ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 4943a40ed3dSBarry Smith PetscFunctionReturn(0); 4958965ea79SLois Curfman McInnes } 4968965ea79SLois Curfman McInnes 4974a2ae208SSatish Balay #undef __FUNCT__ 4984a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 499dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 5008965ea79SLois Curfman McInnes { 5013501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 502dfbe8321SBarry Smith PetscErrorCode ierr; 503ed3cc1f0SBarry Smith 5043a40ed3dSBarry Smith PetscFunctionBegin; 50594d884c6SBarry Smith 506aa482453SBarry Smith #if defined(PETSC_USE_LOG) 50777431f27SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N); 5088965ea79SLois Curfman McInnes #endif 5098798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 510606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 5113501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 5123501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 5133501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 514622d7880SLois Curfman McInnes if (mdn->factor) { 515606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 516606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 517606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 518606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 519622d7880SLois Curfman McInnes } 520606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 521901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr); 522901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 5233a40ed3dSBarry Smith PetscFunctionReturn(0); 5248965ea79SLois Curfman McInnes } 52539ddd567SLois Curfman McInnes 5264a2ae208SSatish Balay #undef __FUNCT__ 5274a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 5286849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5298965ea79SLois Curfman McInnes { 53039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 531dfbe8321SBarry Smith PetscErrorCode ierr; 5327056b6fcSBarry Smith 5333a40ed3dSBarry Smith PetscFunctionBegin; 53439ddd567SLois Curfman McInnes if (mdn->size == 1) { 53539ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5368965ea79SLois Curfman McInnes } 53729bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5383a40ed3dSBarry Smith PetscFunctionReturn(0); 5398965ea79SLois Curfman McInnes } 5408965ea79SLois Curfman McInnes 5414a2ae208SSatish Balay #undef __FUNCT__ 5424a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 5436849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5448965ea79SLois Curfman McInnes { 54539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 546dfbe8321SBarry Smith PetscErrorCode ierr; 54732dcc486SBarry Smith PetscMPIInt size = mdn->size,rank = mdn->rank; 548b0a32e0cSBarry Smith PetscViewerType vtype; 54932077d6dSBarry Smith PetscTruth iascii,isdraw; 550b0a32e0cSBarry Smith PetscViewer sviewer; 551f3ef73ceSBarry Smith PetscViewerFormat format; 5528965ea79SLois Curfman McInnes 5533a40ed3dSBarry Smith PetscFunctionBegin; 55432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 555fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 55632077d6dSBarry Smith if (iascii) { 557b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 558b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 559456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 5604e220ebcSLois Curfman McInnes MatInfo info; 561888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 56277431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m, 56377431f27SBarry Smith (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 564b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5653501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5663a40ed3dSBarry Smith PetscFunctionReturn(0); 567fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 5683a40ed3dSBarry Smith PetscFunctionReturn(0); 5698965ea79SLois Curfman McInnes } 570f1af5d2fSBarry Smith } else if (isdraw) { 571b0a32e0cSBarry Smith PetscDraw draw; 572f1af5d2fSBarry Smith PetscTruth isnull; 573f1af5d2fSBarry Smith 574b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 575b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 576f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 577f1af5d2fSBarry Smith } 57877ed5343SBarry Smith 5798965ea79SLois Curfman McInnes if (size == 1) { 58039ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5813a40ed3dSBarry Smith } else { 5828965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5838965ea79SLois Curfman McInnes Mat A; 58413f74950SBarry Smith PetscInt M = mat->M,N = mat->N,m,row,i,nz; 585ba8c8a56SBarry Smith PetscInt *cols; 586ba8c8a56SBarry Smith PetscScalar *vals; 5878965ea79SLois Curfman McInnes 588*f69a0ea3SMatthew Knepley ierr = MatCreate(mat->comm,&A);CHKERRQ(ierr); 5898965ea79SLois Curfman McInnes if (!rank) { 590*f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 5913a40ed3dSBarry Smith } else { 592*f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 5938965ea79SLois Curfman McInnes } 594be5d1d56SKris Buschelman /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */ 595878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 596878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 59752e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr); 5988965ea79SLois Curfman McInnes 59939ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 60039ddd567SLois Curfman McInnes but it's quick for now */ 60151022da4SBarry Smith A->insertmode = INSERT_VALUES; 602273d9f13SBarry Smith row = mdn->rstart; m = mdn->A->m; 6038965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 604ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 605ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 606ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 60739ddd567SLois Curfman McInnes row++; 6088965ea79SLois Curfman McInnes } 6098965ea79SLois Curfman McInnes 6106d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6116d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 612b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 613b9b97703SBarry Smith if (!rank) { 6146831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 6158965ea79SLois Curfman McInnes } 616b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 617b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6188965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 6198965ea79SLois Curfman McInnes } 6203a40ed3dSBarry Smith PetscFunctionReturn(0); 6218965ea79SLois Curfman McInnes } 6228965ea79SLois Curfman McInnes 6234a2ae208SSatish Balay #undef __FUNCT__ 6244a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 625dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 6268965ea79SLois Curfman McInnes { 627dfbe8321SBarry Smith PetscErrorCode ierr; 62832077d6dSBarry Smith PetscTruth iascii,isbinary,isdraw,issocket; 6298965ea79SLois Curfman McInnes 630433994e6SBarry Smith PetscFunctionBegin; 6310f5bd95cSBarry Smith 63232077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 633fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 634b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 635fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 6360f5bd95cSBarry Smith 63732077d6dSBarry Smith if (iascii || issocket || isdraw) { 638f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 6390f5bd95cSBarry Smith } else if (isbinary) { 6403a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6415cd90555SBarry Smith } else { 642958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6438965ea79SLois Curfman McInnes } 6443a40ed3dSBarry Smith PetscFunctionReturn(0); 6458965ea79SLois Curfman McInnes } 6468965ea79SLois Curfman McInnes 6474a2ae208SSatish Balay #undef __FUNCT__ 6484a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 649dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6508965ea79SLois Curfman McInnes { 6513501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 6523501a2bdSLois Curfman McInnes Mat mdn = mat->A; 653dfbe8321SBarry Smith PetscErrorCode ierr; 654329f5518SBarry Smith PetscReal isend[5],irecv[5]; 6558965ea79SLois Curfman McInnes 6563a40ed3dSBarry Smith PetscFunctionBegin; 657273d9f13SBarry Smith info->rows_global = (double)A->M; 658273d9f13SBarry Smith info->columns_global = (double)A->N; 659273d9f13SBarry Smith info->rows_local = (double)A->m; 660273d9f13SBarry Smith info->columns_local = (double)A->N; 6614e220ebcSLois Curfman McInnes info->block_size = 1.0; 6624e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6634e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6644e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6658965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6664e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6674e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6684e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6694e220ebcSLois Curfman McInnes info->memory = isend[3]; 6704e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6718965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 672d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 6734e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6744e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6754e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6764e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6774e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6788965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 679d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 6804e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6814e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6824e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6834e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6844e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6858965ea79SLois Curfman McInnes } 6864e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6874e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6884e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6893a40ed3dSBarry Smith PetscFunctionReturn(0); 6908965ea79SLois Curfman McInnes } 6918965ea79SLois Curfman McInnes 6924a2ae208SSatish Balay #undef __FUNCT__ 6934a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 694dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op) 6958965ea79SLois Curfman McInnes { 69639ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 697dfbe8321SBarry Smith PetscErrorCode ierr; 6988965ea79SLois Curfman McInnes 6993a40ed3dSBarry Smith PetscFunctionBegin; 70012c028f9SKris Buschelman switch (op) { 70112c028f9SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 70212c028f9SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 70312c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 70412c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 70512c028f9SKris Buschelman case MAT_COLUMNS_SORTED: 70612c028f9SKris Buschelman case MAT_COLUMNS_UNSORTED: 707273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 70812c028f9SKris Buschelman break; 70912c028f9SKris Buschelman case MAT_ROW_ORIENTED: 710273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 711273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 71212c028f9SKris Buschelman break; 71312c028f9SKris Buschelman case MAT_ROWS_SORTED: 71412c028f9SKris Buschelman case MAT_ROWS_UNSORTED: 71512c028f9SKris Buschelman case MAT_YES_NEW_DIAGONALS: 71612c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 71763ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatSetOption_MPIDense:Option ignored\n"));CHKERRQ(ierr); 71812c028f9SKris Buschelman break; 71912c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 720273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 721273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 72212c028f9SKris Buschelman break; 72312c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 724273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 72512c028f9SKris Buschelman break; 72612c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 72729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 72877e54ba9SKris Buschelman case MAT_SYMMETRIC: 72977e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 7309a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 7319a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 7329a4540c5SBarry Smith case MAT_HERMITIAN: 7339a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 7349a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 7359a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 73677e54ba9SKris Buschelman break; 73712c028f9SKris Buschelman default: 73829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 7393a40ed3dSBarry Smith } 7403a40ed3dSBarry Smith PetscFunctionReturn(0); 7418965ea79SLois Curfman McInnes } 7428965ea79SLois Curfman McInnes 7438965ea79SLois Curfman McInnes 7444a2ae208SSatish Balay #undef __FUNCT__ 7454a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 746dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7475b2fa520SLois Curfman McInnes { 7485b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7495b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 75087828ca2SBarry Smith PetscScalar *l,*r,x,*v; 751dfbe8321SBarry Smith PetscErrorCode ierr; 75213f74950SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 7535b2fa520SLois Curfman McInnes 7545b2fa520SLois Curfman McInnes PetscFunctionBegin; 75572d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7565b2fa520SLois Curfman McInnes if (ll) { 75772d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 75829bbc08cSBarry Smith if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 7591ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7605b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7615b2fa520SLois Curfman McInnes x = l[i]; 7625b2fa520SLois Curfman McInnes v = mat->v + i; 7635b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7645b2fa520SLois Curfman McInnes } 7651ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 766efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 7675b2fa520SLois Curfman McInnes } 7685b2fa520SLois Curfman McInnes if (rr) { 76972d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 77029bbc08cSBarry Smith if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 7715b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7725b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7731ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 7745b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 7755b2fa520SLois Curfman McInnes x = r[i]; 7765b2fa520SLois Curfman McInnes v = mat->v + i*m; 7775b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 7785b2fa520SLois Curfman McInnes } 7791ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 780efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 7815b2fa520SLois Curfman McInnes } 7825b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 7835b2fa520SLois Curfman McInnes } 7845b2fa520SLois Curfman McInnes 7854a2ae208SSatish Balay #undef __FUNCT__ 7864a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 787dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 788096963f5SLois Curfman McInnes { 7893501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7903501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 791dfbe8321SBarry Smith PetscErrorCode ierr; 79213f74950SBarry Smith PetscInt i,j; 793329f5518SBarry Smith PetscReal sum = 0.0; 79487828ca2SBarry Smith PetscScalar *v = mat->v; 7953501a2bdSLois Curfman McInnes 7963a40ed3dSBarry Smith PetscFunctionBegin; 7973501a2bdSLois Curfman McInnes if (mdn->size == 1) { 798064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 7993501a2bdSLois Curfman McInnes } else { 8003501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 801273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 802aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 803329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 8043501a2bdSLois Curfman McInnes #else 8053501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 8063501a2bdSLois Curfman McInnes #endif 8073501a2bdSLois Curfman McInnes } 808064f8208SBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 809064f8208SBarry Smith *nrm = sqrt(*nrm); 810efee365bSSatish Balay ierr = PetscLogFlops(2*mdn->A->n*mdn->A->m);CHKERRQ(ierr); 8113a40ed3dSBarry Smith } else if (type == NORM_1) { 812329f5518SBarry Smith PetscReal *tmp,*tmp2; 813b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 814273d9f13SBarry Smith tmp2 = tmp + A->N; 815273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 816064f8208SBarry Smith *nrm = 0.0; 8173501a2bdSLois Curfman McInnes v = mat->v; 818273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 819273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 82067e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8213501a2bdSLois Curfman McInnes } 8223501a2bdSLois Curfman McInnes } 823d7d1e502SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 824273d9f13SBarry Smith for (j=0; j<A->N; j++) { 825064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8263501a2bdSLois Curfman McInnes } 827606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 828efee365bSSatish Balay ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr); 8293a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 830329f5518SBarry Smith PetscReal ntemp; 8313501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 832064f8208SBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 8333a40ed3dSBarry Smith } else { 83429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8353501a2bdSLois Curfman McInnes } 8363501a2bdSLois Curfman McInnes } 8373a40ed3dSBarry Smith PetscFunctionReturn(0); 8383501a2bdSLois Curfman McInnes } 8393501a2bdSLois Curfman McInnes 8404a2ae208SSatish Balay #undef __FUNCT__ 8414a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 842dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout) 8433501a2bdSLois Curfman McInnes { 8443501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8453501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8463501a2bdSLois Curfman McInnes Mat B; 84713f74950SBarry Smith PetscInt M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 8486849ba73SBarry Smith PetscErrorCode ierr; 84913f74950SBarry Smith PetscInt j,i; 85087828ca2SBarry Smith PetscScalar *v; 8513501a2bdSLois Curfman McInnes 8523a40ed3dSBarry Smith PetscFunctionBegin; 8537c922b88SBarry Smith if (!matout && M != N) { 85429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 8557056b6fcSBarry Smith } 856*f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 857*f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr); 858878740d9SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 859878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr); 8603501a2bdSLois Curfman McInnes 861273d9f13SBarry Smith m = a->A->m; n = a->A->n; v = Aloc->v; 86213f74950SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr); 8633501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 8643501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8653501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8663501a2bdSLois Curfman McInnes v += m; 8673501a2bdSLois Curfman McInnes } 868606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 8696d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8706d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8717c922b88SBarry Smith if (matout) { 8723501a2bdSLois Curfman McInnes *matout = B; 8733501a2bdSLois Curfman McInnes } else { 874273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 8753501a2bdSLois Curfman McInnes } 8763a40ed3dSBarry Smith PetscFunctionReturn(0); 877096963f5SLois Curfman McInnes } 878096963f5SLois Curfman McInnes 879d9eff348SSatish Balay #include "petscblaslapack.h" 8804a2ae208SSatish Balay #undef __FUNCT__ 8814a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 882dfbe8321SBarry Smith PetscErrorCode MatScale_MPIDense(const PetscScalar *alpha,Mat inA) 88344cd7ae7SLois Curfman McInnes { 88444cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 88544cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 8864ce68768SBarry Smith PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N; 887efee365bSSatish Balay PetscErrorCode ierr; 88844cd7ae7SLois Curfman McInnes 8893a40ed3dSBarry Smith PetscFunctionBegin; 89071044d3cSBarry Smith BLASscal_(&nz,(PetscScalar*)alpha,a->v,&one); 891efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 8923a40ed3dSBarry Smith PetscFunctionReturn(0); 89344cd7ae7SLois Curfman McInnes } 89444cd7ae7SLois Curfman McInnes 8956849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 8968965ea79SLois Curfman McInnes 8974a2ae208SSatish Balay #undef __FUNCT__ 8984a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 899dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A) 900273d9f13SBarry Smith { 901dfbe8321SBarry Smith PetscErrorCode ierr; 902273d9f13SBarry Smith 903273d9f13SBarry Smith PetscFunctionBegin; 904273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 905273d9f13SBarry Smith PetscFunctionReturn(0); 906273d9f13SBarry Smith } 907273d9f13SBarry Smith 9088965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 90909dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 91009dc0095SBarry Smith MatGetRow_MPIDense, 91109dc0095SBarry Smith MatRestoreRow_MPIDense, 91209dc0095SBarry Smith MatMult_MPIDense, 91397304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 9147c922b88SBarry Smith MatMultTranspose_MPIDense, 9157c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 9168965ea79SLois Curfman McInnes 0, 91709dc0095SBarry Smith 0, 91809dc0095SBarry Smith 0, 91997304618SKris Buschelman /*10*/ 0, 92009dc0095SBarry Smith 0, 92109dc0095SBarry Smith 0, 92209dc0095SBarry Smith 0, 92309dc0095SBarry Smith MatTranspose_MPIDense, 92497304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 92597304618SKris Buschelman 0, 92609dc0095SBarry Smith MatGetDiagonal_MPIDense, 9275b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 92809dc0095SBarry Smith MatNorm_MPIDense, 92997304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 93009dc0095SBarry Smith MatAssemblyEnd_MPIDense, 93109dc0095SBarry Smith 0, 93209dc0095SBarry Smith MatSetOption_MPIDense, 93309dc0095SBarry Smith MatZeroEntries_MPIDense, 93497304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 93509dc0095SBarry Smith 0, 93609dc0095SBarry Smith 0, 93709dc0095SBarry Smith 0, 93809dc0095SBarry Smith 0, 93997304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 940273d9f13SBarry Smith 0, 94109dc0095SBarry Smith 0, 94209dc0095SBarry Smith MatGetArray_MPIDense, 94309dc0095SBarry Smith MatRestoreArray_MPIDense, 94497304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 94509dc0095SBarry Smith 0, 94609dc0095SBarry Smith 0, 94709dc0095SBarry Smith 0, 94809dc0095SBarry Smith 0, 94997304618SKris Buschelman /*40*/ 0, 9502ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 95109dc0095SBarry Smith 0, 95209dc0095SBarry Smith MatGetValues_MPIDense, 95309dc0095SBarry Smith 0, 95497304618SKris Buschelman /*45*/ 0, 95509dc0095SBarry Smith MatScale_MPIDense, 95609dc0095SBarry Smith 0, 95709dc0095SBarry Smith 0, 95809dc0095SBarry Smith 0, 959521d7252SBarry Smith /*50*/ 0, 96009dc0095SBarry Smith 0, 96109dc0095SBarry Smith 0, 96209dc0095SBarry Smith 0, 96309dc0095SBarry Smith 0, 96497304618SKris Buschelman /*55*/ 0, 96509dc0095SBarry Smith 0, 96609dc0095SBarry Smith 0, 96709dc0095SBarry Smith 0, 96809dc0095SBarry Smith 0, 96997304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 970b9b97703SBarry Smith MatDestroy_MPIDense, 971b9b97703SBarry Smith MatView_MPIDense, 97297304618SKris Buschelman MatGetPetscMaps_Petsc, 97397304618SKris Buschelman 0, 97497304618SKris Buschelman /*65*/ 0, 97597304618SKris Buschelman 0, 97697304618SKris Buschelman 0, 97797304618SKris Buschelman 0, 97897304618SKris Buschelman 0, 97997304618SKris Buschelman /*70*/ 0, 98097304618SKris Buschelman 0, 98197304618SKris Buschelman 0, 98297304618SKris Buschelman 0, 98397304618SKris Buschelman 0, 98497304618SKris Buschelman /*75*/ 0, 98597304618SKris Buschelman 0, 98697304618SKris Buschelman 0, 98797304618SKris Buschelman 0, 98897304618SKris Buschelman 0, 98997304618SKris Buschelman /*80*/ 0, 99097304618SKris Buschelman 0, 99197304618SKris Buschelman 0, 99297304618SKris Buschelman 0, 993865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense, 994865e5f61SKris Buschelman 0, 995865e5f61SKris Buschelman 0, 996865e5f61SKris Buschelman 0, 997865e5f61SKris Buschelman 0, 998865e5f61SKris Buschelman 0, 999865e5f61SKris Buschelman /*90*/ 0, 1000865e5f61SKris Buschelman 0, 1001865e5f61SKris Buschelman 0, 1002865e5f61SKris Buschelman 0, 1003865e5f61SKris Buschelman 0, 1004865e5f61SKris Buschelman /*95*/ 0, 1005865e5f61SKris Buschelman 0, 1006865e5f61SKris Buschelman 0, 1007865e5f61SKris Buschelman 0}; 10088965ea79SLois Curfman McInnes 1009273d9f13SBarry Smith EXTERN_C_BEGIN 10104a2ae208SSatish Balay #undef __FUNCT__ 1011a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 1012be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1013a23d5eceSKris Buschelman { 1014a23d5eceSKris Buschelman Mat_MPIDense *a; 1015dfbe8321SBarry Smith PetscErrorCode ierr; 1016a23d5eceSKris Buschelman 1017a23d5eceSKris Buschelman PetscFunctionBegin; 1018a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1019a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1020a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1021a23d5eceSKris Buschelman 1022a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 1023*f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1024*f69a0ea3SMatthew Knepley ierr = MatSetSizes(a->A,mat->m,mat->N,mat->m,mat->N);CHKERRQ(ierr); 10255c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 10265c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 102752e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 1028a23d5eceSKris Buschelman PetscFunctionReturn(0); 1029a23d5eceSKris Buschelman } 1030a23d5eceSKris Buschelman EXTERN_C_END 1031a23d5eceSKris Buschelman 10320bad9183SKris Buschelman /*MC 1033fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 10340bad9183SKris Buschelman 10350bad9183SKris Buschelman Options Database Keys: 10360bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 10370bad9183SKris Buschelman 10380bad9183SKris Buschelman Level: beginner 10390bad9183SKris Buschelman 10400bad9183SKris Buschelman .seealso: MatCreateMPIDense 10410bad9183SKris Buschelman M*/ 10420bad9183SKris Buschelman 1043a23d5eceSKris Buschelman EXTERN_C_BEGIN 1044a23d5eceSKris Buschelman #undef __FUNCT__ 10454a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1046be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat) 1047273d9f13SBarry Smith { 1048273d9f13SBarry Smith Mat_MPIDense *a; 1049dfbe8321SBarry Smith PetscErrorCode ierr; 105013f74950SBarry Smith PetscInt i; 1051273d9f13SBarry Smith 1052273d9f13SBarry Smith PetscFunctionBegin; 1053b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1054b0a32e0cSBarry Smith mat->data = (void*)a; 1055273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1056273d9f13SBarry Smith mat->factor = 0; 1057273d9f13SBarry Smith mat->mapping = 0; 1058273d9f13SBarry Smith 1059273d9f13SBarry Smith a->factor = 0; 1060273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1061273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1062273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1063273d9f13SBarry Smith 1064273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1065273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1066273d9f13SBarry Smith a->nvec = mat->n; 1067273d9f13SBarry Smith 1068273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 1069273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 10708a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 10718a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1072273d9f13SBarry Smith 1073273d9f13SBarry Smith /* build local table of row and column ownerships */ 107413f74950SBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 1075273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 107652e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr); 107713f74950SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1078273d9f13SBarry Smith a->rowners[0] = 0; 1079273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1080273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1081273d9f13SBarry Smith } 1082273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1083273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 108413f74950SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1085273d9f13SBarry Smith a->cowners[0] = 0; 1086273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1087273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1088273d9f13SBarry Smith } 1089273d9f13SBarry Smith 1090273d9f13SBarry Smith /* build cache for off array entries formed */ 1091273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1092273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1093273d9f13SBarry Smith 1094273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1095273d9f13SBarry Smith a->lvec = 0; 1096273d9f13SBarry Smith a->Mvctx = 0; 1097273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1098273d9f13SBarry Smith 1099273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1100273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1101273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1102a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1103a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1104a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1105273d9f13SBarry Smith PetscFunctionReturn(0); 1106273d9f13SBarry Smith } 1107273d9f13SBarry Smith EXTERN_C_END 1108273d9f13SBarry Smith 1109209238afSKris Buschelman /*MC 1110002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1111209238afSKris Buschelman 1112209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1113209238afSKris Buschelman and MATMPIDENSE otherwise. 1114209238afSKris Buschelman 1115209238afSKris Buschelman Options Database Keys: 1116209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1117209238afSKris Buschelman 1118209238afSKris Buschelman Level: beginner 1119209238afSKris Buschelman 1120209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1121209238afSKris Buschelman M*/ 1122209238afSKris Buschelman 1123209238afSKris Buschelman EXTERN_C_BEGIN 1124209238afSKris Buschelman #undef __FUNCT__ 1125209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1126be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A) 1127dfbe8321SBarry Smith { 11286849ba73SBarry Smith PetscErrorCode ierr; 112913f74950SBarry Smith PetscMPIInt size; 1130209238afSKris Buschelman 1131209238afSKris Buschelman PetscFunctionBegin; 1132209238afSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1133209238afSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1134209238afSKris Buschelman if (size == 1) { 1135209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1136209238afSKris Buschelman } else { 1137209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1138209238afSKris Buschelman } 1139209238afSKris Buschelman PetscFunctionReturn(0); 1140209238afSKris Buschelman } 1141209238afSKris Buschelman EXTERN_C_END 1142209238afSKris Buschelman 11434a2ae208SSatish Balay #undef __FUNCT__ 11444a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1145273d9f13SBarry Smith /*@C 1146273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1147273d9f13SBarry Smith 1148273d9f13SBarry Smith Not collective 1149273d9f13SBarry Smith 1150273d9f13SBarry Smith Input Parameters: 1151273d9f13SBarry Smith . A - the matrix 1152273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1153273d9f13SBarry Smith to control all matrix memory allocation. 1154273d9f13SBarry Smith 1155273d9f13SBarry Smith Notes: 1156273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1157273d9f13SBarry Smith storage by columns. 1158273d9f13SBarry Smith 1159273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1160273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1161273d9f13SBarry Smith set data=PETSC_NULL. 1162273d9f13SBarry Smith 1163273d9f13SBarry Smith Level: intermediate 1164273d9f13SBarry Smith 1165273d9f13SBarry Smith .keywords: matrix,dense, parallel 1166273d9f13SBarry Smith 1167273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1168273d9f13SBarry Smith @*/ 1169be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1170273d9f13SBarry Smith { 1171dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar *); 1172273d9f13SBarry Smith 1173273d9f13SBarry Smith PetscFunctionBegin; 1174565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1175a23d5eceSKris Buschelman if (f) { 1176a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1177a23d5eceSKris Buschelman } 1178273d9f13SBarry Smith PetscFunctionReturn(0); 1179273d9f13SBarry Smith } 1180273d9f13SBarry Smith 11814a2ae208SSatish Balay #undef __FUNCT__ 11824a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 11838965ea79SLois Curfman McInnes /*@C 118439ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 11858965ea79SLois Curfman McInnes 1186db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1187db81eaa0SLois Curfman McInnes 11888965ea79SLois Curfman McInnes Input Parameters: 1189db81eaa0SLois Curfman McInnes + comm - MPI communicator 11908965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1191db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 11928965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1193db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 11947f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1195dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 11968965ea79SLois Curfman McInnes 11978965ea79SLois Curfman McInnes Output Parameter: 1198477f1c0bSLois Curfman McInnes . A - the matrix 11998965ea79SLois Curfman McInnes 1200b259b22eSLois Curfman McInnes Notes: 120139ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 120239ddd567SLois Curfman McInnes storage by columns. 12038965ea79SLois Curfman McInnes 120418f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 120518f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 12067f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 120718f449edSLois Curfman McInnes 12088965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 12098965ea79SLois Curfman McInnes (possibly both). 12108965ea79SLois Curfman McInnes 1211027ccd11SLois Curfman McInnes Level: intermediate 1212027ccd11SLois Curfman McInnes 121339ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 12148965ea79SLois Curfman McInnes 121539ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 12168965ea79SLois Curfman McInnes @*/ 1217be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 12188965ea79SLois Curfman McInnes { 12196849ba73SBarry Smith PetscErrorCode ierr; 122013f74950SBarry Smith PetscMPIInt size; 12218965ea79SLois Curfman McInnes 12223a40ed3dSBarry Smith PetscFunctionBegin; 1223*f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1224*f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1225273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1226273d9f13SBarry Smith if (size > 1) { 1227273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1228273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1229273d9f13SBarry Smith } else { 1230273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1231273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 12328c469469SLois Curfman McInnes } 12333a40ed3dSBarry Smith PetscFunctionReturn(0); 12348965ea79SLois Curfman McInnes } 12358965ea79SLois Curfman McInnes 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 12386849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 12398965ea79SLois Curfman McInnes { 12408965ea79SLois Curfman McInnes Mat mat; 12413501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1242dfbe8321SBarry Smith PetscErrorCode ierr; 12438965ea79SLois Curfman McInnes 12443a40ed3dSBarry Smith PetscFunctionBegin; 12458965ea79SLois Curfman McInnes *newmat = 0; 1246*f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&mat);CHKERRQ(ierr); 1247*f69a0ea3SMatthew Knepley ierr = MatSetSizes(mat,A->m,A->n,A->M,A->N);CHKERRQ(ierr); 1248be5d1d56SKris Buschelman ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr); 1249b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1250b0a32e0cSBarry Smith mat->data = (void*)a; 1251549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 12523501a2bdSLois Curfman McInnes mat->factor = A->factor; 1253c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1254273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 12558965ea79SLois Curfman McInnes 12568965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 12578965ea79SLois Curfman McInnes a->rend = oldmat->rend; 12588965ea79SLois Curfman McInnes a->size = oldmat->size; 12598965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1260e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1261b9b97703SBarry Smith a->nvec = oldmat->nvec; 12623782ba37SSatish Balay a->donotstash = oldmat->donotstash; 126313f74950SBarry Smith ierr = PetscMalloc((a->size+1)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 126452e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr); 126513f74950SBarry Smith ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr); 12668798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 12678965ea79SLois Curfman McInnes 1268329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 12695609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 127052e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 12718965ea79SLois Curfman McInnes *newmat = mat; 12723a40ed3dSBarry Smith PetscFunctionReturn(0); 12738965ea79SLois Curfman McInnes } 12748965ea79SLois Curfman McInnes 1275e090d566SSatish Balay #include "petscsys.h" 12768965ea79SLois Curfman McInnes 12774a2ae208SSatish Balay #undef __FUNCT__ 12784a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 1279*f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat) 128090ace30eSBarry Smith { 12816849ba73SBarry Smith PetscErrorCode ierr; 128232dcc486SBarry Smith PetscMPIInt rank,size; 128313f74950SBarry Smith PetscInt *rowners,i,m,nz,j; 128487828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 128590ace30eSBarry Smith MPI_Status status; 128690ace30eSBarry Smith 12873a40ed3dSBarry Smith PetscFunctionBegin; 1288d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1289d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 129090ace30eSBarry Smith 129190ace30eSBarry Smith /* determine ownership of all rows */ 129290ace30eSBarry Smith m = M/size + ((M % size) > rank); 129313f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 129413f74950SBarry Smith ierr = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 129590ace30eSBarry Smith rowners[0] = 0; 129690ace30eSBarry Smith for (i=2; i<=size; i++) { 129790ace30eSBarry Smith rowners[i] += rowners[i-1]; 129890ace30eSBarry Smith } 129990ace30eSBarry Smith 1300*f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1301*f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1302be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1303878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 130490ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 130590ace30eSBarry Smith 130690ace30eSBarry Smith if (!rank) { 130787828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 130890ace30eSBarry Smith 130990ace30eSBarry Smith /* read in my part of the matrix numerical values */ 13100752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 131190ace30eSBarry Smith 131290ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 131390ace30eSBarry Smith vals_ptr = vals; 131490ace30eSBarry Smith for (i=0; i<m; i++) { 131590ace30eSBarry Smith for (j=0; j<N; j++) { 131690ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 131790ace30eSBarry Smith } 131890ace30eSBarry Smith } 131990ace30eSBarry Smith 132090ace30eSBarry Smith /* read in other processors and ship out */ 132190ace30eSBarry Smith for (i=1; i<size; i++) { 132290ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 13230752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1324ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 132590ace30eSBarry Smith } 13263a40ed3dSBarry Smith } else { 132790ace30eSBarry Smith /* receive numeric values */ 132887828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 132990ace30eSBarry Smith 133090ace30eSBarry Smith /* receive message of values*/ 1331ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 133290ace30eSBarry Smith 133390ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 133490ace30eSBarry Smith vals_ptr = vals; 133590ace30eSBarry Smith for (i=0; i<m; i++) { 133690ace30eSBarry Smith for (j=0; j<N; j++) { 133790ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 133890ace30eSBarry Smith } 133990ace30eSBarry Smith } 134090ace30eSBarry Smith } 1341606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1342606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 13436d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13446d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13453a40ed3dSBarry Smith PetscFunctionReturn(0); 134690ace30eSBarry Smith } 134790ace30eSBarry Smith 13484a2ae208SSatish Balay #undef __FUNCT__ 13494a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1350*f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat) 13518965ea79SLois Curfman McInnes { 13528965ea79SLois Curfman McInnes Mat A; 135387828ca2SBarry Smith PetscScalar *vals,*svals; 135419bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 13558965ea79SLois Curfman McInnes MPI_Status status; 135613f74950SBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz; 135713f74950SBarry Smith PetscInt header[4],*rowlengths = 0,M,N,*cols; 135813f74950SBarry Smith PetscInt *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols; 135913f74950SBarry Smith PetscInt i,nz,j,rstart,rend; 136013f74950SBarry Smith int fd; 13616849ba73SBarry Smith PetscErrorCode ierr; 13628965ea79SLois Curfman McInnes 13633a40ed3dSBarry Smith PetscFunctionBegin; 1364d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1365d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13668965ea79SLois Curfman McInnes if (!rank) { 1367b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 13680752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1369552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 13708965ea79SLois Curfman McInnes } 13718965ea79SLois Curfman McInnes 137213f74950SBarry Smith ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 137390ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 137490ace30eSBarry Smith 137590ace30eSBarry Smith /* 137690ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 137790ace30eSBarry Smith */ 137890ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1379be5d1d56SKris Buschelman ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 13803a40ed3dSBarry Smith PetscFunctionReturn(0); 138190ace30eSBarry Smith } 138290ace30eSBarry Smith 13838965ea79SLois Curfman McInnes /* determine ownership of all rows */ 13848965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 138513f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 1386ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 13878965ea79SLois Curfman McInnes rowners[0] = 0; 13888965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 13898965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 13908965ea79SLois Curfman McInnes } 13918965ea79SLois Curfman McInnes rstart = rowners[rank]; 13928965ea79SLois Curfman McInnes rend = rowners[rank+1]; 13938965ea79SLois Curfman McInnes 13948965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 139513f74950SBarry Smith ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr); 13968965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 13978965ea79SLois Curfman McInnes if (!rank) { 139813f74950SBarry Smith ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 13990752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 140013f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr); 14018965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 14021466f1e1SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 1403606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1404ca161407SBarry Smith } else { 14051466f1e1SBarry Smith ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 14068965ea79SLois Curfman McInnes } 14078965ea79SLois Curfman McInnes 14088965ea79SLois Curfman McInnes if (!rank) { 14098965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 141013f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr); 141113f74950SBarry Smith ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 14128965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14138965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 14148965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 14158965ea79SLois Curfman McInnes } 14168965ea79SLois Curfman McInnes } 1417606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 14188965ea79SLois Curfman McInnes 14198965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 14208965ea79SLois Curfman McInnes maxnz = 0; 14218965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14220452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 14238965ea79SLois Curfman McInnes } 142413f74950SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr); 14258965ea79SLois Curfman McInnes 14268965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 14278965ea79SLois Curfman McInnes nz = procsnz[0]; 142813f74950SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14290752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 14308965ea79SLois Curfman McInnes 14318965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 14328965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14338965ea79SLois Curfman McInnes nz = procsnz[i]; 14340752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 143513f74950SBarry Smith ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 14368965ea79SLois Curfman McInnes } 1437606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 14383a40ed3dSBarry Smith } else { 14398965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 14408965ea79SLois Curfman McInnes nz = 0; 14418965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14428965ea79SLois Curfman McInnes nz += ourlens[i]; 14438965ea79SLois Curfman McInnes } 144413f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14458965ea79SLois Curfman McInnes 14468965ea79SLois Curfman McInnes /* receive message of column indices*/ 144713f74950SBarry Smith ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 144813f74950SBarry Smith ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 144929bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14508965ea79SLois Curfman McInnes } 14518965ea79SLois Curfman McInnes 14528965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 145313f74950SBarry Smith ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr); 14548965ea79SLois Curfman McInnes jj = 0; 14558965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14568965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 14578965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 14588965ea79SLois Curfman McInnes jj++; 14598965ea79SLois Curfman McInnes } 14608965ea79SLois Curfman McInnes } 14618965ea79SLois Curfman McInnes 14628965ea79SLois Curfman McInnes /* create our matrix */ 14638965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14648965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 14658965ea79SLois Curfman McInnes } 1466*f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1467*f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1468be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1469878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 14708965ea79SLois Curfman McInnes A = *newmat; 14718965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14728965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 14738965ea79SLois Curfman McInnes } 14748965ea79SLois Curfman McInnes 14758965ea79SLois Curfman McInnes if (!rank) { 147687828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14778965ea79SLois Curfman McInnes 14788965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 14798965ea79SLois Curfman McInnes nz = procsnz[0]; 14800752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 14818965ea79SLois Curfman McInnes 14828965ea79SLois Curfman McInnes /* insert into matrix */ 14838965ea79SLois Curfman McInnes jj = rstart; 14848965ea79SLois Curfman McInnes smycols = mycols; 14858965ea79SLois Curfman McInnes svals = vals; 14868965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14878965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14888965ea79SLois Curfman McInnes smycols += ourlens[i]; 14898965ea79SLois Curfman McInnes svals += ourlens[i]; 14908965ea79SLois Curfman McInnes jj++; 14918965ea79SLois Curfman McInnes } 14928965ea79SLois Curfman McInnes 14938965ea79SLois Curfman McInnes /* read in other processors and ship out */ 14948965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14958965ea79SLois Curfman McInnes nz = procsnz[i]; 14960752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1497ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 14988965ea79SLois Curfman McInnes } 1499606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 15003a40ed3dSBarry Smith } else { 15018965ea79SLois Curfman McInnes /* receive numeric values */ 150284d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 15038965ea79SLois Curfman McInnes 15048965ea79SLois Curfman McInnes /* receive message of values*/ 1505ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1506ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 150729bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 15088965ea79SLois Curfman McInnes 15098965ea79SLois Curfman McInnes /* insert into matrix */ 15108965ea79SLois Curfman McInnes jj = rstart; 15118965ea79SLois Curfman McInnes smycols = mycols; 15128965ea79SLois Curfman McInnes svals = vals; 15138965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15148965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 15158965ea79SLois Curfman McInnes smycols += ourlens[i]; 15168965ea79SLois Curfman McInnes svals += ourlens[i]; 15178965ea79SLois Curfman McInnes jj++; 15188965ea79SLois Curfman McInnes } 15198965ea79SLois Curfman McInnes } 1520606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1521606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1522606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1523606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 15248965ea79SLois Curfman McInnes 15256d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15266d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15273a40ed3dSBarry Smith PetscFunctionReturn(0); 15288965ea79SLois Curfman McInnes } 152990ace30eSBarry Smith 153090ace30eSBarry Smith 153190ace30eSBarry Smith 153290ace30eSBarry Smith 153390ace30eSBarry Smith 1534