1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3ed3cc1f0SBarry Smith /* 4ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 5ed3cc1f0SBarry Smith */ 6ed3cc1f0SBarry Smith 704fea9ffSBarry Smith 8*7c4f633dSBarry Smith #include "../src/mat/impls/dense/mpi/mpidense.h" /*I "petscmat.h" I*/ 9eae6fb2eSBarry Smith #if defined(PETSC_HAVE_PLAPACK) 10eae6fb2eSBarry Smith static PetscMPIInt Plapack_nprows,Plapack_npcols,Plapack_ierror,Plapack_nb_alg; 11eae6fb2eSBarry Smith static MPI_Comm Plapack_comm_2d; 12eae6fb2eSBarry Smith #endif 138965ea79SLois Curfman McInnes 14ba8c8a56SBarry Smith #undef __FUNCT__ 15ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix" 16ab92ecdeSBarry Smith /*@ 17ab92ecdeSBarry Smith 18ab92ecdeSBarry Smith MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential 19ab92ecdeSBarry Smith matrix that represents the operator. For sequential matrices it returns itself. 20ab92ecdeSBarry Smith 21ab92ecdeSBarry Smith Input Parameter: 22ab92ecdeSBarry Smith . A - the Seq or MPI dense matrix 23ab92ecdeSBarry Smith 24ab92ecdeSBarry Smith Output Parameter: 25ab92ecdeSBarry Smith . B - the inner matrix 26ab92ecdeSBarry Smith 278e6c10adSSatish Balay Level: intermediate 288e6c10adSSatish Balay 29ab92ecdeSBarry Smith @*/ 30ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B) 31ab92ecdeSBarry Smith { 32ab92ecdeSBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 33ab92ecdeSBarry Smith PetscErrorCode ierr; 34ab92ecdeSBarry Smith PetscTruth flg; 35ab92ecdeSBarry Smith 36ab92ecdeSBarry Smith PetscFunctionBegin; 37ab92ecdeSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr); 38ab92ecdeSBarry Smith if (flg) { 39ab92ecdeSBarry Smith *B = mat->A; 40ab92ecdeSBarry Smith } else { 41ab92ecdeSBarry Smith *B = A; 42ab92ecdeSBarry Smith } 43ab92ecdeSBarry Smith PetscFunctionReturn(0); 44ab92ecdeSBarry Smith } 45ab92ecdeSBarry Smith 46ab92ecdeSBarry Smith #undef __FUNCT__ 47ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense" 48ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 49ba8c8a56SBarry Smith { 50ba8c8a56SBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 51ba8c8a56SBarry Smith PetscErrorCode ierr; 52d0f46423SBarry Smith PetscInt lrow,rstart = A->rmap->rstart,rend = A->rmap->rend; 53ba8c8a56SBarry Smith 54ba8c8a56SBarry Smith PetscFunctionBegin; 55ba8c8a56SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 56ba8c8a56SBarry Smith lrow = row - rstart; 57ba8c8a56SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr); 58ba8c8a56SBarry Smith PetscFunctionReturn(0); 59ba8c8a56SBarry Smith } 60ba8c8a56SBarry Smith 61ba8c8a56SBarry Smith #undef __FUNCT__ 62ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense" 63ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 64ba8c8a56SBarry Smith { 65ba8c8a56SBarry Smith PetscErrorCode ierr; 66ba8c8a56SBarry Smith 67ba8c8a56SBarry Smith PetscFunctionBegin; 68ba8c8a56SBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 69ba8c8a56SBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 70ba8c8a56SBarry Smith PetscFunctionReturn(0); 71ba8c8a56SBarry Smith } 72ba8c8a56SBarry Smith 730de54da6SSatish Balay EXTERN_C_BEGIN 744a2ae208SSatish Balay #undef __FUNCT__ 754a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 76be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 770de54da6SSatish Balay { 780de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 796849ba73SBarry Smith PetscErrorCode ierr; 80d0f46423SBarry Smith PetscInt m = A->rmap->n,rstart = A->rmap->rstart; 8187828ca2SBarry Smith PetscScalar *array; 820de54da6SSatish Balay MPI_Comm comm; 830de54da6SSatish Balay 840de54da6SSatish Balay PetscFunctionBegin; 85d0f46423SBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 860de54da6SSatish Balay 870de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 880de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 890de54da6SSatish Balay 900de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 910de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 92f69a0ea3SMatthew Knepley ierr = MatCreate(comm,B);CHKERRQ(ierr); 93f69a0ea3SMatthew Knepley ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr); 947adad957SLisandro Dalcin ierr = MatSetType(*B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr); 955c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr); 960de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 970de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 980de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 990de54da6SSatish Balay 1000de54da6SSatish Balay *iscopy = PETSC_TRUE; 1010de54da6SSatish Balay PetscFunctionReturn(0); 1020de54da6SSatish Balay } 1030de54da6SSatish Balay EXTERN_C_END 1040de54da6SSatish Balay 1054a2ae208SSatish Balay #undef __FUNCT__ 1064a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 10713f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) 1088965ea79SLois Curfman McInnes { 10939b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 110dfbe8321SBarry Smith PetscErrorCode ierr; 111d0f46423SBarry Smith PetscInt i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row; 112273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 1138965ea79SLois Curfman McInnes 1143a40ed3dSBarry Smith PetscFunctionBegin; 1158965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 1165ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 117d0f46423SBarry Smith if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 1188965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 1198965ea79SLois Curfman McInnes row = idxm[i] - rstart; 12039b7565bSBarry Smith if (roworiented) { 12139b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 1223a40ed3dSBarry Smith } else { 1238965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 1245ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 125d0f46423SBarry Smith if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 12639b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 12739b7565bSBarry Smith } 1288965ea79SLois Curfman McInnes } 1293a40ed3dSBarry Smith } else { 1303782ba37SSatish Balay if (!A->donotstash) { 13139b7565bSBarry Smith if (roworiented) { 1328798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 133d36fbae8SSatish Balay } else { 1348798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 13539b7565bSBarry Smith } 136b49de8d1SLois Curfman McInnes } 137b49de8d1SLois Curfman McInnes } 1383782ba37SSatish Balay } 1393a40ed3dSBarry Smith PetscFunctionReturn(0); 140b49de8d1SLois Curfman McInnes } 141b49de8d1SLois Curfman McInnes 1424a2ae208SSatish Balay #undef __FUNCT__ 1434a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 14413f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 145b49de8d1SLois Curfman McInnes { 146b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 147dfbe8321SBarry Smith PetscErrorCode ierr; 148d0f46423SBarry Smith PetscInt i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row; 149b49de8d1SLois Curfman McInnes 1503a40ed3dSBarry Smith PetscFunctionBegin; 151b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 15297e567efSBarry Smith if (idxm[i] < 0) continue; /* SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */ 153d0f46423SBarry Smith if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 154b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 155b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 156b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 15797e567efSBarry Smith if (idxn[j] < 0) continue; /* SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */ 158d0f46423SBarry Smith if (idxn[j] >= mat->cmap->N) { 15929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 160a8c6a408SBarry Smith } 161b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 162b49de8d1SLois Curfman McInnes } 163a8c6a408SBarry Smith } else { 16429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 1658965ea79SLois Curfman McInnes } 1668965ea79SLois Curfman McInnes } 1673a40ed3dSBarry Smith PetscFunctionReturn(0); 1688965ea79SLois Curfman McInnes } 1698965ea79SLois Curfman McInnes 1704a2ae208SSatish Balay #undef __FUNCT__ 1714a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 172dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 173ff14e315SSatish Balay { 174ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 175dfbe8321SBarry Smith PetscErrorCode ierr; 176ff14e315SSatish Balay 1773a40ed3dSBarry Smith PetscFunctionBegin; 178ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1793a40ed3dSBarry Smith PetscFunctionReturn(0); 180ff14e315SSatish Balay } 181ff14e315SSatish Balay 1824a2ae208SSatish Balay #undef __FUNCT__ 1834a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 18413f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 185ca3fa75bSLois Curfman McInnes { 186ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 187ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 1886849ba73SBarry Smith PetscErrorCode ierr; 1895d0c19d7SBarry Smith PetscInt i,j,rstart,rend,nrows,ncols,nlrows,nlcols; 1905d0c19d7SBarry Smith const PetscInt *irow,*icol; 19187828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 192ca3fa75bSLois Curfman McInnes Mat newmat; 193ca3fa75bSLois Curfman McInnes 194ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 195ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 196ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 197b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 198b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 199ca3fa75bSLois Curfman McInnes 200ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 2017eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 2027eba5e9cSLois Curfman McInnes original matrix! */ 203ca3fa75bSLois Curfman McInnes 204ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 2057eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 206ca3fa75bSLois Curfman McInnes 207ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 208ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 20929bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 2107eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 211ca3fa75bSLois Curfman McInnes newmat = *B; 212ca3fa75bSLois Curfman McInnes } else { 213ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 2147adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 215f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr); 2167adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 217878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 218ca3fa75bSLois Curfman McInnes } 219ca3fa75bSLois Curfman McInnes 220ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 221ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 222ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 223ca3fa75bSLois Curfman McInnes 224ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 225ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 226ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 2277eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 228ca3fa75bSLois Curfman McInnes } 229ca3fa75bSLois Curfman McInnes } 230ca3fa75bSLois Curfman McInnes 231ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 232ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 233ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 234ca3fa75bSLois Curfman McInnes 235ca3fa75bSLois Curfman McInnes /* Free work space */ 236ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 237ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 238ca3fa75bSLois Curfman McInnes *B = newmat; 239ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 240ca3fa75bSLois Curfman McInnes } 241ca3fa75bSLois Curfman McInnes 2424a2ae208SSatish Balay #undef __FUNCT__ 2434a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 244dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 245ff14e315SSatish Balay { 2463a40ed3dSBarry Smith PetscFunctionBegin; 2473a40ed3dSBarry Smith PetscFunctionReturn(0); 248ff14e315SSatish Balay } 249ff14e315SSatish Balay 2504a2ae208SSatish Balay #undef __FUNCT__ 2514a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 252dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2538965ea79SLois Curfman McInnes { 25439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 2557adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)mat)->comm; 256dfbe8321SBarry Smith PetscErrorCode ierr; 25713f74950SBarry Smith PetscInt nstash,reallocs; 2588965ea79SLois Curfman McInnes InsertMode addv; 2598965ea79SLois Curfman McInnes 2603a40ed3dSBarry Smith PetscFunctionBegin; 2618965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 262ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 2637056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 26429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 2658965ea79SLois Curfman McInnes } 266e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2678965ea79SLois Curfman McInnes 268d0f46423SBarry Smith ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr); 2698798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 270ae15b995SBarry Smith ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr); 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 2728965ea79SLois Curfman McInnes } 2738965ea79SLois Curfman McInnes 2744a2ae208SSatish Balay #undef __FUNCT__ 2754a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 276dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2778965ea79SLois Curfman McInnes { 27839ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2796849ba73SBarry Smith PetscErrorCode ierr; 28013f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 28113f74950SBarry Smith PetscMPIInt n; 28287828ca2SBarry Smith PetscScalar *val; 283e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2848965ea79SLois Curfman McInnes 2853a40ed3dSBarry Smith PetscFunctionBegin; 2868965ea79SLois Curfman McInnes /* wait on receives */ 2877ef1d9bdSSatish Balay while (1) { 2888798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2897ef1d9bdSSatish Balay if (!flg) break; 2908965ea79SLois Curfman McInnes 2917ef1d9bdSSatish Balay for (i=0; i<n;) { 2927ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2937ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2947ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2957ef1d9bdSSatish Balay else ncols = n-i; 2967ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2977ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2987ef1d9bdSSatish Balay i = j; 2998965ea79SLois Curfman McInnes } 3007ef1d9bdSSatish Balay } 3018798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 3028965ea79SLois Curfman McInnes 30339ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 30439ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 3058965ea79SLois Curfman McInnes 3066d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 30739ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 3088965ea79SLois Curfman McInnes } 3093a40ed3dSBarry Smith PetscFunctionReturn(0); 3108965ea79SLois Curfman McInnes } 3118965ea79SLois Curfman McInnes 3124a2ae208SSatish Balay #undef __FUNCT__ 3134a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 314dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 3158965ea79SLois Curfman McInnes { 316dfbe8321SBarry Smith PetscErrorCode ierr; 31739ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3183a40ed3dSBarry Smith 3193a40ed3dSBarry Smith PetscFunctionBegin; 3203a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 3213a40ed3dSBarry Smith PetscFunctionReturn(0); 3228965ea79SLois Curfman McInnes } 3238965ea79SLois Curfman McInnes 3248965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 3258965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 3268965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 3273501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 3288965ea79SLois Curfman McInnes routine. 3298965ea79SLois Curfman McInnes */ 3304a2ae208SSatish Balay #undef __FUNCT__ 3314a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 332f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 3338965ea79SLois Curfman McInnes { 33439ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3356849ba73SBarry Smith PetscErrorCode ierr; 336d0f46423SBarry Smith PetscInt i,*owners = A->rmap->range; 33713f74950SBarry Smith PetscInt *nprocs,j,idx,nsends; 33813f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 3397adad957SLisandro Dalcin PetscInt *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source; 34013f74950SBarry Smith PetscInt *lens,*lrows,*values; 34113f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 3427adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)A)->comm; 3438965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3448965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 34535d8aa7fSBarry Smith PetscTruth found; 3468965ea79SLois Curfman McInnes 3473a40ed3dSBarry Smith PetscFunctionBegin; 3488965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 34913f74950SBarry Smith ierr = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr); 35013f74950SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 35113f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/ 3528965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3538965ea79SLois Curfman McInnes idx = rows[i]; 35435d8aa7fSBarry Smith found = PETSC_FALSE; 3558965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3568965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 357c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3588965ea79SLois Curfman McInnes } 3598965ea79SLois Curfman McInnes } 36029bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3618965ea79SLois Curfman McInnes } 362c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 3638965ea79SLois Curfman McInnes 3648965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 365c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3668965ea79SLois Curfman McInnes 3678965ea79SLois Curfman McInnes /* post receives: */ 36813f74950SBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr); 369b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3708965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 37113f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3728965ea79SLois Curfman McInnes } 3738965ea79SLois Curfman McInnes 3748965ea79SLois Curfman McInnes /* do sends: 3758965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3768965ea79SLois Curfman McInnes the ith processor 3778965ea79SLois Curfman McInnes */ 37813f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr); 379b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 38013f74950SBarry Smith ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr); 3818965ea79SLois Curfman McInnes starts[0] = 0; 382c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3838965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3848965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3858965ea79SLois Curfman McInnes } 3868965ea79SLois Curfman McInnes 3878965ea79SLois Curfman McInnes starts[0] = 0; 388c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3898965ea79SLois Curfman McInnes count = 0; 3908965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 391c1dc657dSBarry Smith if (nprocs[2*i+1]) { 39213f74950SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3938965ea79SLois Curfman McInnes } 3948965ea79SLois Curfman McInnes } 395606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3968965ea79SLois Curfman McInnes 3978965ea79SLois Curfman McInnes base = owners[rank]; 3988965ea79SLois Curfman McInnes 3998965ea79SLois Curfman McInnes /* wait on receives */ 40013f74950SBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr); 4018965ea79SLois Curfman McInnes source = lens + nrecvs; 4028965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 4038965ea79SLois Curfman McInnes while (count) { 404ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4058965ea79SLois Curfman McInnes /* unpack receives into our local space */ 40613f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 4078965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 4088965ea79SLois Curfman McInnes lens[imdex] = n; 4098965ea79SLois Curfman McInnes slen += n; 4108965ea79SLois Curfman McInnes count--; 4118965ea79SLois Curfman McInnes } 412606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 4138965ea79SLois Curfman McInnes 4148965ea79SLois Curfman McInnes /* move the data into the send scatter */ 41513f74950SBarry Smith ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr); 4168965ea79SLois Curfman McInnes count = 0; 4178965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 4188965ea79SLois Curfman McInnes values = rvalues + i*nmax; 4198965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 4208965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 4218965ea79SLois Curfman McInnes } 4228965ea79SLois Curfman McInnes } 423606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 424606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 425606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 426606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 4278965ea79SLois Curfman McInnes 4288965ea79SLois Curfman McInnes /* actually zap the local rows */ 429f4df32b1SMatthew Knepley ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr); 430606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 4318965ea79SLois Curfman McInnes 4328965ea79SLois Curfman McInnes /* wait on sends */ 4338965ea79SLois Curfman McInnes if (nsends) { 434b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 435ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 436606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 4378965ea79SLois Curfman McInnes } 438606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 439606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 4408965ea79SLois Curfman McInnes 4413a40ed3dSBarry Smith PetscFunctionReturn(0); 4428965ea79SLois Curfman McInnes } 4438965ea79SLois Curfman McInnes 4444a2ae208SSatish Balay #undef __FUNCT__ 4454a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 446dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4478965ea79SLois Curfman McInnes { 44839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 449dfbe8321SBarry Smith PetscErrorCode ierr; 450c456f294SBarry Smith 4513a40ed3dSBarry Smith PetscFunctionBegin; 452ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 453ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 45444cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4553a40ed3dSBarry Smith PetscFunctionReturn(0); 4568965ea79SLois Curfman McInnes } 4578965ea79SLois Curfman McInnes 4584a2ae208SSatish Balay #undef __FUNCT__ 4594a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 460dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4618965ea79SLois Curfman McInnes { 46239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 463dfbe8321SBarry Smith PetscErrorCode ierr; 464c456f294SBarry Smith 4653a40ed3dSBarry Smith PetscFunctionBegin; 466ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 467ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 46844cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4693a40ed3dSBarry Smith PetscFunctionReturn(0); 4708965ea79SLois Curfman McInnes } 4718965ea79SLois Curfman McInnes 4724a2ae208SSatish Balay #undef __FUNCT__ 4734a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 474dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 475096963f5SLois Curfman McInnes { 476096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 477dfbe8321SBarry Smith PetscErrorCode ierr; 47887828ca2SBarry Smith PetscScalar zero = 0.0; 479096963f5SLois Curfman McInnes 4803a40ed3dSBarry Smith PetscFunctionBegin; 4812dcb1b2aSMatthew Knepley ierr = VecSet(yy,zero);CHKERRQ(ierr); 4827c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 483ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 484ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4853a40ed3dSBarry Smith PetscFunctionReturn(0); 486096963f5SLois Curfman McInnes } 487096963f5SLois Curfman McInnes 4884a2ae208SSatish Balay #undef __FUNCT__ 4894a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 490dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 491096963f5SLois Curfman McInnes { 492096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 493dfbe8321SBarry Smith PetscErrorCode ierr; 494096963f5SLois Curfman McInnes 4953a40ed3dSBarry Smith PetscFunctionBegin; 4963501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4977c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 498ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 499ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5003a40ed3dSBarry Smith PetscFunctionReturn(0); 501096963f5SLois Curfman McInnes } 502096963f5SLois Curfman McInnes 5034a2ae208SSatish Balay #undef __FUNCT__ 5044a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 505dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 5068965ea79SLois Curfman McInnes { 50739ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 508096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 509dfbe8321SBarry Smith PetscErrorCode ierr; 510d0f46423SBarry Smith PetscInt len,i,n,m = A->rmap->n,radd; 51187828ca2SBarry Smith PetscScalar *x,zero = 0.0; 512ed3cc1f0SBarry Smith 5133a40ed3dSBarry Smith PetscFunctionBegin; 5142dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 5151ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 516096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 517d0f46423SBarry Smith if (n != A->rmap->N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 518d0f46423SBarry Smith len = PetscMin(a->A->rmap->n,a->A->cmap->n); 519d0f46423SBarry Smith radd = A->rmap->rstart*m; 52044cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 521096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 522096963f5SLois Curfman McInnes } 5231ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 5243a40ed3dSBarry Smith PetscFunctionReturn(0); 5258965ea79SLois Curfman McInnes } 5268965ea79SLois Curfman McInnes 5274a2ae208SSatish Balay #undef __FUNCT__ 5284a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 529dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 5308965ea79SLois Curfman McInnes { 5313501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 532dfbe8321SBarry Smith PetscErrorCode ierr; 53301b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK) 53401b82886SBarry Smith Mat_Plapack *lu=(Mat_Plapack*)(mat->spptr); 53501b82886SBarry Smith #endif 536ed3cc1f0SBarry Smith 5373a40ed3dSBarry Smith PetscFunctionBegin; 53894d884c6SBarry Smith 539aa482453SBarry Smith #if defined(PETSC_USE_LOG) 540d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N); 5418965ea79SLois Curfman McInnes #endif 5428798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 5433501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 54405b42c5fSBarry Smith if (mdn->lvec) {ierr = VecDestroy(mdn->lvec);CHKERRQ(ierr);} 54505b42c5fSBarry Smith if (mdn->Mvctx) {ierr = VecScatterDestroy(mdn->Mvctx);CHKERRQ(ierr);} 54601b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK) 5472fbe02b9SBarry Smith if (lu) { 54862b4c0b3SBarry Smith ierr = PLA_Obj_free(&lu->A);CHKERRQ(ierr); 54962b4c0b3SBarry Smith ierr = PLA_Obj_free (&lu->pivots);CHKERRQ(ierr); 55062b4c0b3SBarry Smith ierr = PLA_Temp_free(&lu->templ);CHKERRQ(ierr); 55109d27a7eSBarry Smith 5522fbe02b9SBarry Smith if (lu->is_pla) { 55301b82886SBarry Smith ierr = ISDestroy(lu->is_pla);CHKERRQ(ierr); 55401b82886SBarry Smith ierr = ISDestroy(lu->is_petsc);CHKERRQ(ierr); 55501b82886SBarry Smith ierr = VecScatterDestroy(lu->ctx);CHKERRQ(ierr); 556622d7880SLois Curfman McInnes } 5572fbe02b9SBarry Smith } 55801b82886SBarry Smith #endif 55901b82886SBarry Smith 560606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 561dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 562901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr); 563901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 5644ae313f4SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr); 5654ae313f4SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr); 5664ae313f4SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr); 5673a40ed3dSBarry Smith PetscFunctionReturn(0); 5688965ea79SLois Curfman McInnes } 56939ddd567SLois Curfman McInnes 5704a2ae208SSatish Balay #undef __FUNCT__ 5714a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 5726849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5738965ea79SLois Curfman McInnes { 57439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 575dfbe8321SBarry Smith PetscErrorCode ierr; 576aa05aa95SBarry Smith PetscViewerFormat format; 577aa05aa95SBarry Smith int fd; 578d0f46423SBarry Smith PetscInt header[4],mmax,N = mat->cmap->N,i,j,m,k; 579aa05aa95SBarry Smith PetscMPIInt rank,tag = ((PetscObject)viewer)->tag,size; 580578230a0SSatish Balay PetscScalar *work,*v,*vv; 581aa05aa95SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)mdn->A->data; 582aa05aa95SBarry Smith MPI_Status status; 5837056b6fcSBarry Smith 5843a40ed3dSBarry Smith PetscFunctionBegin; 58539ddd567SLois Curfman McInnes if (mdn->size == 1) { 58639ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 587aa05aa95SBarry Smith } else { 588aa05aa95SBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 5897adad957SLisandro Dalcin ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&rank);CHKERRQ(ierr); 5907adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr); 591aa05aa95SBarry Smith 592aa05aa95SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 5937973ad04SBarry Smith if (format == PETSC_VIEWER_NATIVE) { 594aa05aa95SBarry Smith 595aa05aa95SBarry Smith if (!rank) { 596aa05aa95SBarry Smith /* store the matrix as a dense matrix */ 597aa05aa95SBarry Smith header[0] = MAT_FILE_COOKIE; 598d0f46423SBarry Smith header[1] = mat->rmap->N; 599aa05aa95SBarry Smith header[2] = N; 600aa05aa95SBarry Smith header[3] = MATRIX_BINARY_FORMAT_DENSE; 601aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 602aa05aa95SBarry Smith 603aa05aa95SBarry Smith /* get largest work array needed for transposing array */ 604d0f46423SBarry Smith mmax = mat->rmap->n; 605aa05aa95SBarry Smith for (i=1; i<size; i++) { 606d0f46423SBarry Smith mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]); 6078965ea79SLois Curfman McInnes } 608aa05aa95SBarry Smith ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&work);CHKERRQ(ierr); 609aa05aa95SBarry Smith 610aa05aa95SBarry Smith /* write out local array, by rows */ 611d0f46423SBarry Smith m = mat->rmap->n; 612aa05aa95SBarry Smith v = a->v; 613aa05aa95SBarry Smith for (j=0; j<N; j++) { 614aa05aa95SBarry Smith for (i=0; i<m; i++) { 615578230a0SSatish Balay work[j + i*N] = *v++; 616aa05aa95SBarry Smith } 617aa05aa95SBarry Smith } 618aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 619aa05aa95SBarry Smith /* get largest work array to receive messages from other processes, excludes process zero */ 620aa05aa95SBarry Smith mmax = 0; 621aa05aa95SBarry Smith for (i=1; i<size; i++) { 622d0f46423SBarry Smith mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]); 623aa05aa95SBarry Smith } 624578230a0SSatish Balay ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&vv);CHKERRQ(ierr); 625aa05aa95SBarry Smith for(k = 1; k < size; k++) { 626f8009846SMatthew Knepley v = vv; 627d0f46423SBarry Smith m = mat->rmap->range[k+1] - mat->rmap->range[k]; 6287adad957SLisandro Dalcin ierr = MPI_Recv(v,m*N,MPIU_SCALAR,k,tag,((PetscObject)mat)->comm,&status);CHKERRQ(ierr); 629aa05aa95SBarry Smith 630aa05aa95SBarry Smith for(j = 0; j < N; j++) { 631aa05aa95SBarry Smith for(i = 0; i < m; i++) { 632578230a0SSatish Balay work[j + i*N] = *v++; 633aa05aa95SBarry Smith } 634aa05aa95SBarry Smith } 635aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 636aa05aa95SBarry Smith } 637aa05aa95SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 638578230a0SSatish Balay ierr = PetscFree(vv);CHKERRQ(ierr); 639aa05aa95SBarry Smith } else { 640d0f46423SBarry Smith ierr = MPI_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,((PetscObject)mat)->comm);CHKERRQ(ierr); 641aa05aa95SBarry Smith } 6424aa34b0aSBarry Smith } else { 6437973ad04SBarry Smith SETERRQ(PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE"); 644aa05aa95SBarry Smith } 645aa05aa95SBarry Smith } 6463a40ed3dSBarry Smith PetscFunctionReturn(0); 6478965ea79SLois Curfman McInnes } 6488965ea79SLois Curfman McInnes 6494a2ae208SSatish Balay #undef __FUNCT__ 6504a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 6516849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 6528965ea79SLois Curfman McInnes { 65339ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 654dfbe8321SBarry Smith PetscErrorCode ierr; 65532dcc486SBarry Smith PetscMPIInt size = mdn->size,rank = mdn->rank; 656a313700dSBarry Smith const PetscViewerType vtype; 65732077d6dSBarry Smith PetscTruth iascii,isdraw; 658b0a32e0cSBarry Smith PetscViewer sviewer; 659f3ef73ceSBarry Smith PetscViewerFormat format; 66001b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK) 66101b82886SBarry Smith Mat_Plapack *lu=(Mat_Plapack*)(mat->spptr); 66201b82886SBarry Smith #endif 6638965ea79SLois Curfman McInnes 6643a40ed3dSBarry Smith PetscFunctionBegin; 66532077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 666fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 66732077d6dSBarry Smith if (iascii) { 668b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 669b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 670456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 6714e220ebcSLois Curfman McInnes MatInfo info; 672888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 673d0f46423SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap->n, 67477431f27SBarry Smith (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 675b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6763501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 67701b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK) 67801b82886SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PLAPACK run parameters:\n");CHKERRQ(ierr); 6796c995c7dSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Processor mesh: nprows %d, npcols %d\n",lu->Plapack_nprows, lu->Plapack_npcols);CHKERRQ(ierr); 68001b82886SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Distr. block size nb: %d \n",lu->nb);CHKERRQ(ierr); 6816c995c7dSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Error checking: %d\n",lu->Plapack_ierror);CHKERRQ(ierr); 6826c995c7dSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Algorithmic block size: %d\n",lu->Plapack_nb_alg);CHKERRQ(ierr); 68301b82886SBarry Smith #endif 6843a40ed3dSBarry Smith PetscFunctionReturn(0); 685fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 6863a40ed3dSBarry Smith PetscFunctionReturn(0); 6878965ea79SLois Curfman McInnes } 688f1af5d2fSBarry Smith } else if (isdraw) { 689b0a32e0cSBarry Smith PetscDraw draw; 690f1af5d2fSBarry Smith PetscTruth isnull; 691f1af5d2fSBarry Smith 692b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 693b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 694f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 695f1af5d2fSBarry Smith } 69677ed5343SBarry Smith 6978965ea79SLois Curfman McInnes if (size == 1) { 69839ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 6993a40ed3dSBarry Smith } else { 7008965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 7018965ea79SLois Curfman McInnes Mat A; 702d0f46423SBarry Smith PetscInt M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz; 703ba8c8a56SBarry Smith PetscInt *cols; 704ba8c8a56SBarry Smith PetscScalar *vals; 7058965ea79SLois Curfman McInnes 7067adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)mat)->comm,&A);CHKERRQ(ierr); 7078965ea79SLois Curfman McInnes if (!rank) { 708f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 7093a40ed3dSBarry Smith } else { 710f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 7118965ea79SLois Curfman McInnes } 7127adad957SLisandro Dalcin /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */ 713878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 714878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 71552e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr); 7168965ea79SLois Curfman McInnes 71739ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 71839ddd567SLois Curfman McInnes but it's quick for now */ 71951022da4SBarry Smith A->insertmode = INSERT_VALUES; 720d0f46423SBarry Smith row = mat->rmap->rstart; m = mdn->A->rmap->n; 7218965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 722ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 723ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 724ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 72539ddd567SLois Curfman McInnes row++; 7268965ea79SLois Curfman McInnes } 7278965ea79SLois Curfman McInnes 7286d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7296d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 730b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 731b9b97703SBarry Smith if (!rank) { 7326831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 7338965ea79SLois Curfman McInnes } 734b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 735b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 7368965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 7378965ea79SLois Curfman McInnes } 7383a40ed3dSBarry Smith PetscFunctionReturn(0); 7398965ea79SLois Curfman McInnes } 7408965ea79SLois Curfman McInnes 7414a2ae208SSatish Balay #undef __FUNCT__ 7424a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 743dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 7448965ea79SLois Curfman McInnes { 745dfbe8321SBarry Smith PetscErrorCode ierr; 74632077d6dSBarry Smith PetscTruth iascii,isbinary,isdraw,issocket; 7478965ea79SLois Curfman McInnes 748433994e6SBarry Smith PetscFunctionBegin; 7490f5bd95cSBarry Smith 75032077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 751fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 752b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 753fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 7540f5bd95cSBarry Smith 75532077d6dSBarry Smith if (iascii || issocket || isdraw) { 756f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 7570f5bd95cSBarry Smith } else if (isbinary) { 7583a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 7595cd90555SBarry Smith } else { 760958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 7618965ea79SLois Curfman McInnes } 7623a40ed3dSBarry Smith PetscFunctionReturn(0); 7638965ea79SLois Curfman McInnes } 7648965ea79SLois Curfman McInnes 7654a2ae208SSatish Balay #undef __FUNCT__ 7664a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 767dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 7688965ea79SLois Curfman McInnes { 7693501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7703501a2bdSLois Curfman McInnes Mat mdn = mat->A; 771dfbe8321SBarry Smith PetscErrorCode ierr; 772329f5518SBarry Smith PetscReal isend[5],irecv[5]; 7738965ea79SLois Curfman McInnes 7743a40ed3dSBarry Smith PetscFunctionBegin; 7754e220ebcSLois Curfman McInnes info->block_size = 1.0; 7764e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 7774e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 7784e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 7798965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 7804e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 7814e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 7824e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 7834e220ebcSLois Curfman McInnes info->memory = isend[3]; 7844e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 7858965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 7867adad957SLisandro Dalcin ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,((PetscObject)A)->comm);CHKERRQ(ierr); 7874e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7884e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7894e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7904e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7914e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7928965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 7937adad957SLisandro Dalcin ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 7944e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7954e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7964e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7974e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7984e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7998965ea79SLois Curfman McInnes } 8004e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 8014e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 8024e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 8033a40ed3dSBarry Smith PetscFunctionReturn(0); 8048965ea79SLois Curfman McInnes } 8058965ea79SLois Curfman McInnes 8064a2ae208SSatish Balay #undef __FUNCT__ 8074a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 8084e0d8c25SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscTruth flg) 8098965ea79SLois Curfman McInnes { 81039ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 811dfbe8321SBarry Smith PetscErrorCode ierr; 8128965ea79SLois Curfman McInnes 8133a40ed3dSBarry Smith PetscFunctionBegin; 81412c028f9SKris Buschelman switch (op) { 815512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 81612c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 81712c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 8184e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 81912c028f9SKris Buschelman break; 82012c028f9SKris Buschelman case MAT_ROW_ORIENTED: 8214e0d8c25SBarry Smith a->roworiented = flg; 8224e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 82312c028f9SKris Buschelman break; 8244e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 82512c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 826290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 82712c028f9SKris Buschelman break; 82812c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 8294e0d8c25SBarry Smith a->donotstash = flg; 83012c028f9SKris Buschelman break; 83177e54ba9SKris Buschelman case MAT_SYMMETRIC: 83277e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 8339a4540c5SBarry Smith case MAT_HERMITIAN: 8349a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 835600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 836290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 83777e54ba9SKris Buschelman break; 83812c028f9SKris Buschelman default: 839600fe468SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 8403a40ed3dSBarry Smith } 8413a40ed3dSBarry Smith PetscFunctionReturn(0); 8428965ea79SLois Curfman McInnes } 8438965ea79SLois Curfman McInnes 8448965ea79SLois Curfman McInnes 8454a2ae208SSatish Balay #undef __FUNCT__ 8464a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 847dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 8485b2fa520SLois Curfman McInnes { 8495b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8505b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 85187828ca2SBarry Smith PetscScalar *l,*r,x,*v; 852dfbe8321SBarry Smith PetscErrorCode ierr; 853d0f46423SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n; 8545b2fa520SLois Curfman McInnes 8555b2fa520SLois Curfman McInnes PetscFunctionBegin; 85672d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 8575b2fa520SLois Curfman McInnes if (ll) { 85872d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 859175be7b4SMatthew Knepley if (s2a != s2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2); 8601ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 8615b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 8625b2fa520SLois Curfman McInnes x = l[i]; 8635b2fa520SLois Curfman McInnes v = mat->v + i; 8645b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 8655b2fa520SLois Curfman McInnes } 8661ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 867efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8685b2fa520SLois Curfman McInnes } 8695b2fa520SLois Curfman McInnes if (rr) { 870175be7b4SMatthew Knepley ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr); 871175be7b4SMatthew Knepley if (s3a != s3) SETERRQ2(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3); 872ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 873ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8741ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 8755b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 8765b2fa520SLois Curfman McInnes x = r[i]; 8775b2fa520SLois Curfman McInnes v = mat->v + i*m; 8785b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 8795b2fa520SLois Curfman McInnes } 8801ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 881efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8825b2fa520SLois Curfman McInnes } 8835b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 8845b2fa520SLois Curfman McInnes } 8855b2fa520SLois Curfman McInnes 8864a2ae208SSatish Balay #undef __FUNCT__ 8874a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 888dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 889096963f5SLois Curfman McInnes { 8903501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8913501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 892dfbe8321SBarry Smith PetscErrorCode ierr; 89313f74950SBarry Smith PetscInt i,j; 894329f5518SBarry Smith PetscReal sum = 0.0; 89587828ca2SBarry Smith PetscScalar *v = mat->v; 8963501a2bdSLois Curfman McInnes 8973a40ed3dSBarry Smith PetscFunctionBegin; 8983501a2bdSLois Curfman McInnes if (mdn->size == 1) { 899064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 9003501a2bdSLois Curfman McInnes } else { 9013501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 902d0f46423SBarry Smith for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) { 903aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 904329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 9053501a2bdSLois Curfman McInnes #else 9063501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 9073501a2bdSLois Curfman McInnes #endif 9083501a2bdSLois Curfman McInnes } 9097adad957SLisandro Dalcin ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 910064f8208SBarry Smith *nrm = sqrt(*nrm); 911d0f46423SBarry Smith ierr = PetscLogFlops(2*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr); 9123a40ed3dSBarry Smith } else if (type == NORM_1) { 913329f5518SBarry Smith PetscReal *tmp,*tmp2; 914d0f46423SBarry Smith ierr = PetscMalloc(2*A->cmap->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 915d0f46423SBarry Smith tmp2 = tmp + A->cmap->N; 916d0f46423SBarry Smith ierr = PetscMemzero(tmp,2*A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr); 917064f8208SBarry Smith *nrm = 0.0; 9183501a2bdSLois Curfman McInnes v = mat->v; 919d0f46423SBarry Smith for (j=0; j<mdn->A->cmap->n; j++) { 920d0f46423SBarry Smith for (i=0; i<mdn->A->rmap->n; i++) { 92167e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 9223501a2bdSLois Curfman McInnes } 9233501a2bdSLois Curfman McInnes } 924d0f46423SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr); 925d0f46423SBarry Smith for (j=0; j<A->cmap->N; j++) { 926064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 9273501a2bdSLois Curfman McInnes } 928606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 929d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 9303a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 931329f5518SBarry Smith PetscReal ntemp; 9323501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 9337adad957SLisandro Dalcin ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,((PetscObject)A)->comm);CHKERRQ(ierr); 9343a40ed3dSBarry Smith } else { 93529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 9363501a2bdSLois Curfman McInnes } 9373501a2bdSLois Curfman McInnes } 9383a40ed3dSBarry Smith PetscFunctionReturn(0); 9393501a2bdSLois Curfman McInnes } 9403501a2bdSLois Curfman McInnes 9414a2ae208SSatish Balay #undef __FUNCT__ 9424a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 943fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout) 9443501a2bdSLois Curfman McInnes { 9453501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 9463501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 9473501a2bdSLois Curfman McInnes Mat B; 948d0f46423SBarry Smith PetscInt M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart; 9496849ba73SBarry Smith PetscErrorCode ierr; 95013f74950SBarry Smith PetscInt j,i; 95187828ca2SBarry Smith PetscScalar *v; 9523501a2bdSLois Curfman McInnes 9533a40ed3dSBarry Smith PetscFunctionBegin; 954e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 955fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX || A == *matout) { 9567adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 957d0f46423SBarry Smith ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr); 9587adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 959878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr); 960fc4dec0aSBarry Smith } else { 961fc4dec0aSBarry Smith B = *matout; 962fc4dec0aSBarry Smith } 9633501a2bdSLois Curfman McInnes 964d0f46423SBarry Smith m = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v; 9651acff37aSSatish Balay ierr = PetscMalloc(m*sizeof(PetscInt),&rwork);CHKERRQ(ierr); 9663501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 9671acff37aSSatish Balay for (j=0; j<n; j++) { 9683501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 9693501a2bdSLois Curfman McInnes v += m; 9703501a2bdSLois Curfman McInnes } 971606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 9726d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9736d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 974815cbec1SBarry Smith if (reuse == MAT_INITIAL_MATRIX || *matout != A) { 9753501a2bdSLois Curfman McInnes *matout = B; 9763501a2bdSLois Curfman McInnes } else { 977273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 9783501a2bdSLois Curfman McInnes } 9793a40ed3dSBarry Smith PetscFunctionReturn(0); 980096963f5SLois Curfman McInnes } 981096963f5SLois Curfman McInnes 982d9eff348SSatish Balay #include "petscblaslapack.h" 9834a2ae208SSatish Balay #undef __FUNCT__ 9844a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 985f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha) 98644cd7ae7SLois Curfman McInnes { 98744cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 98844cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 989f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 990efee365bSSatish Balay PetscErrorCode ierr; 991d0f46423SBarry Smith PetscBLASInt one = 1,nz = PetscBLASIntCast(inA->rmap->n*inA->cmap->N); 99244cd7ae7SLois Curfman McInnes 9933a40ed3dSBarry Smith PetscFunctionBegin; 994f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 995efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 9963a40ed3dSBarry Smith PetscFunctionReturn(0); 99744cd7ae7SLois Curfman McInnes } 99844cd7ae7SLois Curfman McInnes 9996849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 10008965ea79SLois Curfman McInnes 10014a2ae208SSatish Balay #undef __FUNCT__ 10024a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 1003dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A) 1004273d9f13SBarry Smith { 1005dfbe8321SBarry Smith PetscErrorCode ierr; 1006273d9f13SBarry Smith 1007273d9f13SBarry Smith PetscFunctionBegin; 1008273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 1009273d9f13SBarry Smith PetscFunctionReturn(0); 1010273d9f13SBarry Smith } 1011273d9f13SBarry Smith 101201b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK) 1013eae6fb2eSBarry Smith 1014eae6fb2eSBarry Smith #undef __FUNCT__ 1015eae6fb2eSBarry Smith #define __FUNCT__ "MatMPIDenseCopyToPlapack" 1016eae6fb2eSBarry Smith PetscErrorCode MatMPIDenseCopyToPlapack(Mat A,Mat F) 1017eae6fb2eSBarry Smith { 1018eae6fb2eSBarry Smith Mat_Plapack *lu = (Mat_Plapack*)(F)->spptr; 1019eae6fb2eSBarry Smith PetscErrorCode ierr; 1020d0f46423SBarry Smith PetscInt M=A->cmap->N,m=A->rmap->n,rstart; 1021eae6fb2eSBarry Smith PetscScalar *array; 1022eae6fb2eSBarry Smith PetscReal one = 1.0; 1023eae6fb2eSBarry Smith 1024eae6fb2eSBarry Smith PetscFunctionBegin; 10252fbe02b9SBarry Smith /* Copy A into F->lu->A */ 1026eae6fb2eSBarry Smith ierr = PLA_Obj_set_to_zero(lu->A);CHKERRQ(ierr); 1027eae6fb2eSBarry Smith ierr = PLA_API_begin();CHKERRQ(ierr); 1028eae6fb2eSBarry Smith ierr = PLA_Obj_API_open(lu->A);CHKERRQ(ierr); 102979b0a62dSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 1030eae6fb2eSBarry Smith ierr = MatGetArray(A,&array);CHKERRQ(ierr); 1031eae6fb2eSBarry Smith ierr = PLA_API_axpy_matrix_to_global(m,M, &one,(void *)array,m,lu->A,rstart,0);CHKERRQ(ierr); 1032eae6fb2eSBarry Smith ierr = MatRestoreArray(A,&array);CHKERRQ(ierr); 1033eae6fb2eSBarry Smith ierr = PLA_Obj_API_close(lu->A);CHKERRQ(ierr); 1034eae6fb2eSBarry Smith ierr = PLA_API_end();CHKERRQ(ierr); 1035eae6fb2eSBarry Smith lu->rstart = rstart; 1036eae6fb2eSBarry Smith PetscFunctionReturn(0); 1037eae6fb2eSBarry Smith } 1038eae6fb2eSBarry Smith 103901b82886SBarry Smith #undef __FUNCT__ 10402fbe02b9SBarry Smith #define __FUNCT__ "MatMPIDenseCopyFromPlapack" 10412fbe02b9SBarry Smith PetscErrorCode MatMPIDenseCopyFromPlapack(Mat F,Mat A) 10422fbe02b9SBarry Smith { 10432fbe02b9SBarry Smith Mat_Plapack *lu = (Mat_Plapack*)(F)->spptr; 10442fbe02b9SBarry Smith PetscErrorCode ierr; 1045d0f46423SBarry Smith PetscInt M=A->cmap->N,m=A->rmap->n,rstart; 10462fbe02b9SBarry Smith PetscScalar *array; 10472fbe02b9SBarry Smith PetscReal one = 1.0; 10482fbe02b9SBarry Smith 10492fbe02b9SBarry Smith PetscFunctionBegin; 10502fbe02b9SBarry Smith /* Copy F into A->lu->A */ 105179b0a62dSBarry Smith ierr = MatZeroEntries(A);CHKERRQ(ierr); 10522fbe02b9SBarry Smith ierr = PLA_API_begin();CHKERRQ(ierr); 10532fbe02b9SBarry Smith ierr = PLA_Obj_API_open(lu->A);CHKERRQ(ierr); 10542fbe02b9SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 10552fbe02b9SBarry Smith ierr = MatGetArray(A,&array);CHKERRQ(ierr); 10562fbe02b9SBarry Smith ierr = PLA_API_axpy_global_to_matrix(m,M, &one,lu->A,rstart,0,(void *)array,m);CHKERRQ(ierr); 10572fbe02b9SBarry Smith ierr = MatRestoreArray(A,&array);CHKERRQ(ierr); 10582fbe02b9SBarry Smith ierr = PLA_Obj_API_close(lu->A);CHKERRQ(ierr); 10592fbe02b9SBarry Smith ierr = PLA_API_end();CHKERRQ(ierr); 10602fbe02b9SBarry Smith lu->rstart = rstart; 10612fbe02b9SBarry Smith PetscFunctionReturn(0); 10622fbe02b9SBarry Smith } 10632fbe02b9SBarry Smith 10642fbe02b9SBarry Smith #undef __FUNCT__ 10652fbe02b9SBarry Smith #define __FUNCT__ "MatMatMultNumeric_MPIDense_MPIDense" 10662fbe02b9SBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C) 10672fbe02b9SBarry Smith { 10682fbe02b9SBarry Smith PetscErrorCode ierr; 10692fbe02b9SBarry Smith Mat_Plapack *luA = (Mat_Plapack*)A->spptr; 10702fbe02b9SBarry Smith Mat_Plapack *luB = (Mat_Plapack*)B->spptr; 10712fbe02b9SBarry Smith Mat_Plapack *luC = (Mat_Plapack*)C->spptr; 10722fbe02b9SBarry Smith PLA_Obj alpha = NULL,beta = NULL; 10732fbe02b9SBarry Smith 10742fbe02b9SBarry Smith PetscFunctionBegin; 10752fbe02b9SBarry Smith ierr = MatMPIDenseCopyToPlapack(A,A);CHKERRQ(ierr); 10762fbe02b9SBarry Smith ierr = MatMPIDenseCopyToPlapack(B,B);CHKERRQ(ierr); 10772fbe02b9SBarry Smith 10786e6f9017SBarry Smith /* 1079cb2480eaSBarry Smith ierr = PLA_Global_show("A = ",luA->A,"%g ","");CHKERRQ(ierr); 1080cb2480eaSBarry Smith ierr = PLA_Global_show("B = ",luB->A,"%g ","");CHKERRQ(ierr); 10816e6f9017SBarry Smith */ 1082cb2480eaSBarry Smith 10832fbe02b9SBarry Smith /* do the multiply in PLA */ 108479b0a62dSBarry Smith ierr = PLA_Create_constants_conf_to(luA->A,NULL,NULL,&alpha);CHKERRQ(ierr); 108579b0a62dSBarry Smith ierr = PLA_Create_constants_conf_to(luC->A,NULL,&beta,NULL);CHKERRQ(ierr); 108679b0a62dSBarry Smith CHKMEMQ; 10872fbe02b9SBarry Smith 1088f0e3846dSSatish Balay ierr = PLA_Gemm(PLA_NO_TRANSPOSE,PLA_NO_TRANSPOSE,alpha,luA->A,luB->A,beta,luC->A); /* CHKERRQ(ierr); */ 108979b0a62dSBarry Smith CHKMEMQ; 10902fbe02b9SBarry Smith ierr = PLA_Obj_free(&alpha);CHKERRQ(ierr); 10912fbe02b9SBarry Smith ierr = PLA_Obj_free(&beta);CHKERRQ(ierr); 10922fbe02b9SBarry Smith 10936e6f9017SBarry Smith /* 1094cb2480eaSBarry Smith ierr = PLA_Global_show("C = ",luC->A,"%g ","");CHKERRQ(ierr); 10956e6f9017SBarry Smith */ 10962fbe02b9SBarry Smith ierr = MatMPIDenseCopyFromPlapack(C,C);CHKERRQ(ierr); 10972fbe02b9SBarry Smith PetscFunctionReturn(0); 10982fbe02b9SBarry Smith } 10992fbe02b9SBarry Smith 11002fbe02b9SBarry Smith #undef __FUNCT__ 11012fbe02b9SBarry Smith #define __FUNCT__ "MatMatMultSymbolic_MPIDense_MPIDense" 11022fbe02b9SBarry Smith PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 11032fbe02b9SBarry Smith { 11042fbe02b9SBarry Smith PetscErrorCode ierr; 1105d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 11062fbe02b9SBarry Smith Mat Cmat; 11072fbe02b9SBarry Smith 11082fbe02b9SBarry Smith PetscFunctionBegin; 1109d0f46423SBarry Smith if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n); 11106e6f9017SBarry Smith SETERRQ(PETSC_ERR_LIB,"Due to aparent bugs in PLAPACK,this is not currently supported"); 11112fbe02b9SBarry Smith ierr = MatCreate(((PetscObject)B)->comm,&Cmat);CHKERRQ(ierr); 1112d0f46423SBarry Smith ierr = MatSetSizes(Cmat,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 11132fbe02b9SBarry Smith ierr = MatSetType(Cmat,MATMPIDENSE);CHKERRQ(ierr); 11142fbe02b9SBarry Smith ierr = MatAssemblyBegin(Cmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11152fbe02b9SBarry Smith ierr = MatAssemblyEnd(Cmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 111685bd4cefSHong Zhang 11172fbe02b9SBarry Smith *C = Cmat; 11182fbe02b9SBarry Smith PetscFunctionReturn(0); 11192fbe02b9SBarry Smith } 11202fbe02b9SBarry Smith 11212fbe02b9SBarry Smith #undef __FUNCT__ 11222fbe02b9SBarry Smith #define __FUNCT__ "MatMatMult_MPIDense_MPIDense" 11232fbe02b9SBarry Smith PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 11242fbe02b9SBarry Smith { 11252fbe02b9SBarry Smith PetscErrorCode ierr; 11262fbe02b9SBarry Smith 11272fbe02b9SBarry Smith PetscFunctionBegin; 11282fbe02b9SBarry Smith if (scall == MAT_INITIAL_MATRIX){ 11292fbe02b9SBarry Smith ierr = MatMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr); 11302fbe02b9SBarry Smith } 11312fbe02b9SBarry Smith ierr = MatMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr); 11322fbe02b9SBarry Smith PetscFunctionReturn(0); 11332fbe02b9SBarry Smith } 11342fbe02b9SBarry Smith 1135ab235cb6SHong Zhang #undef __FUNCT__ 11366c995c7dSHong Zhang #define __FUNCT__ "MatSolve_MPIDense" 11376c995c7dSHong Zhang PetscErrorCode MatSolve_MPIDense(Mat A,Vec b,Vec x) 11386c995c7dSHong Zhang { 11396c995c7dSHong Zhang MPI_Comm comm = ((PetscObject)A)->comm; 11406c995c7dSHong Zhang Mat_Plapack *lu = (Mat_Plapack*)A->spptr; 11416c995c7dSHong Zhang PetscErrorCode ierr; 1142d0f46423SBarry Smith PetscInt M=A->rmap->N,m=A->rmap->n,rstart,i,j,*idx_pla,*idx_petsc,loc_m,loc_stride; 11436c995c7dSHong Zhang PetscScalar *array; 11446c995c7dSHong Zhang PetscReal one = 1.0; 11456c995c7dSHong Zhang PetscMPIInt size,rank,r_rank,r_nproc,c_rank,c_nproc;; 11466c995c7dSHong Zhang PLA_Obj v_pla = NULL; 11476c995c7dSHong Zhang PetscScalar *loc_buf; 11486c995c7dSHong Zhang Vec loc_x; 11496c995c7dSHong Zhang 11506c995c7dSHong Zhang PetscFunctionBegin; 11516c995c7dSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 11526c995c7dSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 11536c995c7dSHong Zhang 11546c995c7dSHong Zhang /* Create PLAPACK vector objects, then copy b into PLAPACK b */ 11556c995c7dSHong Zhang PLA_Mvector_create(lu->datatype,M,1,lu->templ,PLA_ALIGN_FIRST,&v_pla); 11566c995c7dSHong Zhang PLA_Obj_set_to_zero(v_pla); 11576c995c7dSHong Zhang 11586c995c7dSHong Zhang /* Copy b into rhs_pla */ 11596c995c7dSHong Zhang PLA_API_begin(); 11606c995c7dSHong Zhang PLA_Obj_API_open(v_pla); 11616c995c7dSHong Zhang ierr = VecGetArray(b,&array);CHKERRQ(ierr); 11626c995c7dSHong Zhang PLA_API_axpy_vector_to_global(m,&one,(void *)array,1,v_pla,lu->rstart); 11636c995c7dSHong Zhang ierr = VecRestoreArray(b,&array);CHKERRQ(ierr); 11646c995c7dSHong Zhang PLA_Obj_API_close(v_pla); 11656c995c7dSHong Zhang PLA_API_end(); 11666c995c7dSHong Zhang 11676c995c7dSHong Zhang if (A->factor == MAT_FACTOR_LU){ 11686c995c7dSHong Zhang /* Apply the permutations to the right hand sides */ 11696c995c7dSHong Zhang PLA_Apply_pivots_to_rows (v_pla,lu->pivots); 11706c995c7dSHong Zhang 11716c995c7dSHong Zhang /* Solve L y = b, overwriting b with y */ 11726c995c7dSHong Zhang PLA_Trsv( PLA_LOWER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_UNIT_DIAG,lu->A,v_pla ); 11736c995c7dSHong Zhang 11746c995c7dSHong Zhang /* Solve U x = y (=b), overwriting b with x */ 11756c995c7dSHong Zhang PLA_Trsv( PLA_UPPER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_NONUNIT_DIAG,lu->A,v_pla ); 11766c995c7dSHong Zhang } else { /*MAT_FACTOR_CHOLESKY */ 11776c995c7dSHong Zhang PLA_Trsv( PLA_LOWER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_NONUNIT_DIAG,lu->A,v_pla); 11786c995c7dSHong Zhang PLA_Trsv( PLA_LOWER_TRIANGULAR,(lu->datatype == MPI_DOUBLE ? PLA_TRANSPOSE : PLA_CONJUGATE_TRANSPOSE), 11796c995c7dSHong Zhang PLA_NONUNIT_DIAG,lu->A,v_pla); 11806c995c7dSHong Zhang } 11816c995c7dSHong Zhang 11826c995c7dSHong Zhang /* Copy PLAPACK x into Petsc vector x */ 11836c995c7dSHong Zhang PLA_Obj_local_length(v_pla, &loc_m); 11846c995c7dSHong Zhang PLA_Obj_local_buffer(v_pla, (void**)&loc_buf); 11856c995c7dSHong Zhang PLA_Obj_local_stride(v_pla, &loc_stride); 11866c995c7dSHong Zhang /* 11876c995c7dSHong Zhang PetscPrintf(PETSC_COMM_SELF," [%d] b - local_m %d local_stride %d, loc_buf: %g %g, nb: %d\n",rank,loc_m,loc_stride,loc_buf[0],loc_buf[(loc_m-1)*loc_stride],lu->nb); 11886c995c7dSHong Zhang */ 11896c995c7dSHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,loc_m*loc_stride,loc_buf,&loc_x);CHKERRQ(ierr); 11906c995c7dSHong Zhang if (!lu->pla_solved){ 11916c995c7dSHong Zhang 11926c995c7dSHong Zhang PLA_Temp_comm_row_info(lu->templ,&lu->comm_2d,&r_rank,&r_nproc); 11936c995c7dSHong Zhang PLA_Temp_comm_col_info(lu->templ,&lu->comm_2d,&c_rank,&c_nproc); 11946c995c7dSHong Zhang /* printf(" [%d] rank: %d %d, nproc: %d %d\n",rank,r_rank,c_rank,r_nproc,c_nproc); */ 11956c995c7dSHong Zhang 11966c995c7dSHong Zhang /* Create IS and cts for VecScatterring */ 11976c995c7dSHong Zhang PLA_Obj_local_length(v_pla, &loc_m); 11986c995c7dSHong Zhang PLA_Obj_local_stride(v_pla, &loc_stride); 11996c995c7dSHong Zhang ierr = PetscMalloc((2*loc_m+1)*sizeof(PetscInt),&idx_pla);CHKERRQ(ierr); 12006c995c7dSHong Zhang idx_petsc = idx_pla + loc_m; 12016c995c7dSHong Zhang 12026c995c7dSHong Zhang rstart = (r_rank*c_nproc+c_rank)*lu->nb; 12036c995c7dSHong Zhang for (i=0; i<loc_m; i+=lu->nb){ 12046c995c7dSHong Zhang j = 0; 12056c995c7dSHong Zhang while (j < lu->nb && i+j < loc_m){ 12066c995c7dSHong Zhang idx_petsc[i+j] = rstart + j; j++; 12076c995c7dSHong Zhang } 12086c995c7dSHong Zhang rstart += size*lu->nb; 12096c995c7dSHong Zhang } 12106c995c7dSHong Zhang 12116c995c7dSHong Zhang for (i=0; i<loc_m; i++) idx_pla[i] = i*loc_stride; 12126c995c7dSHong Zhang 12136c995c7dSHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,loc_m,idx_pla,&lu->is_pla);CHKERRQ(ierr); 12146c995c7dSHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,loc_m,idx_petsc,&lu->is_petsc);CHKERRQ(ierr); 12156c995c7dSHong Zhang ierr = PetscFree(idx_pla);CHKERRQ(ierr); 12166c995c7dSHong Zhang ierr = VecScatterCreate(loc_x,lu->is_pla,x,lu->is_petsc,&lu->ctx);CHKERRQ(ierr); 12176c995c7dSHong Zhang } 12186c995c7dSHong Zhang ierr = VecScatterBegin(lu->ctx,loc_x,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12196c995c7dSHong Zhang ierr = VecScatterEnd(lu->ctx,loc_x,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12206c995c7dSHong Zhang 12216c995c7dSHong Zhang /* Free data */ 12226c995c7dSHong Zhang ierr = VecDestroy(loc_x);CHKERRQ(ierr); 12236c995c7dSHong Zhang PLA_Obj_free(&v_pla); 12246c995c7dSHong Zhang 12256c995c7dSHong Zhang lu->pla_solved = PETSC_TRUE; 12266c995c7dSHong Zhang PetscFunctionReturn(0); 12276c995c7dSHong Zhang } 12286c995c7dSHong Zhang 12296c995c7dSHong Zhang #undef __FUNCT__ 12306c995c7dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_MPIDense" 12310481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_MPIDense(Mat F,Mat A,const MatFactorInfo *info) 12326c995c7dSHong Zhang { 1233719d5645SBarry Smith Mat_Plapack *lu = (Mat_Plapack*)(F)->spptr; 12346c995c7dSHong Zhang PetscErrorCode ierr; 1235d0f46423SBarry Smith PetscInt M=A->rmap->N,m=A->rmap->n,rstart,rend; 12366c995c7dSHong Zhang PetscInt info_pla=0; 12376c995c7dSHong Zhang PetscScalar *array,one = 1.0; 12386c995c7dSHong Zhang 12396c995c7dSHong Zhang PetscFunctionBegin; 12406c995c7dSHong Zhang if (lu->mstruct == SAME_NONZERO_PATTERN){ 12416c995c7dSHong Zhang PLA_Obj_free(&lu->A); 12426c995c7dSHong Zhang PLA_Obj_free (&lu->pivots); 12436c995c7dSHong Zhang } 12446c995c7dSHong Zhang /* Create PLAPACK matrix object */ 12456c995c7dSHong Zhang lu->A = NULL; lu->pivots = NULL; 12466c995c7dSHong Zhang PLA_Matrix_create(lu->datatype,M,M,lu->templ,PLA_ALIGN_FIRST,PLA_ALIGN_FIRST,&lu->A); 12476c995c7dSHong Zhang PLA_Obj_set_to_zero(lu->A); 12486c995c7dSHong Zhang PLA_Mvector_create(MPI_INT,M,1,lu->templ,PLA_ALIGN_FIRST,&lu->pivots); 12496c995c7dSHong Zhang 12506c995c7dSHong Zhang /* Copy A into lu->A */ 12516c995c7dSHong Zhang PLA_API_begin(); 12526c995c7dSHong Zhang PLA_Obj_API_open(lu->A); 12536c995c7dSHong Zhang ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 12546c995c7dSHong Zhang ierr = MatGetArray(A,&array);CHKERRQ(ierr); 12556c995c7dSHong Zhang PLA_API_axpy_matrix_to_global(m,M, &one,(void *)array,m,lu->A,rstart,0); 12566c995c7dSHong Zhang ierr = MatRestoreArray(A,&array);CHKERRQ(ierr); 12576c995c7dSHong Zhang PLA_Obj_API_close(lu->A); 12586c995c7dSHong Zhang PLA_API_end(); 12596c995c7dSHong Zhang 12606c995c7dSHong Zhang /* Factor P A -> L U overwriting lower triangular portion of A with L, upper, U */ 12616c995c7dSHong Zhang info_pla = PLA_LU(lu->A,lu->pivots); 12626c995c7dSHong Zhang if (info_pla != 0) 12636c995c7dSHong Zhang SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot encountered at row %d from PLA_LU()",info_pla); 12646c995c7dSHong Zhang 12656c995c7dSHong Zhang lu->CleanUpPlapack = PETSC_TRUE; 12666c995c7dSHong Zhang lu->rstart = rstart; 12676c995c7dSHong Zhang lu->mstruct = SAME_NONZERO_PATTERN; 1268719d5645SBarry Smith F->ops->solve = MatSolve_MPIDense; 1269719d5645SBarry Smith F->assembled = PETSC_TRUE; /* required by -ksp_view */ 12706c995c7dSHong Zhang PetscFunctionReturn(0); 12716c995c7dSHong Zhang } 12726c995c7dSHong Zhang 12736c995c7dSHong Zhang #undef __FUNCT__ 12746c995c7dSHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_MPIDense" 12750481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_MPIDense(Mat F,Mat A,const MatFactorInfo *info) 12766c995c7dSHong Zhang { 1277719d5645SBarry Smith Mat_Plapack *lu = (Mat_Plapack*)F->spptr; 12786c995c7dSHong Zhang PetscErrorCode ierr; 1279d0f46423SBarry Smith PetscInt M=A->rmap->N,m=A->rmap->n,rstart,rend; 12806c995c7dSHong Zhang PetscInt info_pla=0; 12816c995c7dSHong Zhang PetscScalar *array,one = 1.0; 12826c995c7dSHong Zhang 12836c995c7dSHong Zhang PetscFunctionBegin; 12846c995c7dSHong Zhang if (lu->mstruct == SAME_NONZERO_PATTERN){ 12856c995c7dSHong Zhang PLA_Obj_free(&lu->A); 12866c995c7dSHong Zhang } 12876c995c7dSHong Zhang /* Create PLAPACK matrix object */ 12886c995c7dSHong Zhang lu->A = NULL; 12896c995c7dSHong Zhang lu->pivots = NULL; 12906c995c7dSHong Zhang PLA_Matrix_create(lu->datatype,M,M,lu->templ,PLA_ALIGN_FIRST,PLA_ALIGN_FIRST,&lu->A); 12916c995c7dSHong Zhang 12926c995c7dSHong Zhang /* Copy A into lu->A */ 12936c995c7dSHong Zhang PLA_API_begin(); 12946c995c7dSHong Zhang PLA_Obj_API_open(lu->A); 12956c995c7dSHong Zhang ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 12966c995c7dSHong Zhang ierr = MatGetArray(A,&array);CHKERRQ(ierr); 12976c995c7dSHong Zhang PLA_API_axpy_matrix_to_global(m,M, &one,(void *)array,m,lu->A,rstart,0); 12986c995c7dSHong Zhang ierr = MatRestoreArray(A,&array);CHKERRQ(ierr); 12996c995c7dSHong Zhang PLA_Obj_API_close(lu->A); 13006c995c7dSHong Zhang PLA_API_end(); 13016c995c7dSHong Zhang 13026c995c7dSHong Zhang /* Factor P A -> Chol */ 13036c995c7dSHong Zhang info_pla = PLA_Chol(PLA_LOWER_TRIANGULAR,lu->A); 13046c995c7dSHong Zhang if (info_pla != 0) 13056c995c7dSHong Zhang SETERRQ1( PETSC_ERR_MAT_CH_ZRPVT,"Nonpositive definite matrix detected at row %d from PLA_Chol()",info_pla); 13066c995c7dSHong Zhang 13076c995c7dSHong Zhang lu->CleanUpPlapack = PETSC_TRUE; 13086c995c7dSHong Zhang lu->rstart = rstart; 13096c995c7dSHong Zhang lu->mstruct = SAME_NONZERO_PATTERN; 1310719d5645SBarry Smith F->ops->solve = MatSolve_MPIDense; 1311719d5645SBarry Smith F->assembled = PETSC_TRUE; /* required by -ksp_view */ 13126c995c7dSHong Zhang PetscFunctionReturn(0); 13136c995c7dSHong Zhang } 13146c995c7dSHong Zhang 13156c995c7dSHong Zhang #undef __FUNCT__ 1316ab235cb6SHong Zhang #define __FUNCT__ "MatFactorSymbolic_Plapack_Private" 13170481f469SBarry Smith PetscErrorCode MatFactorSymbolic_Plapack_Private(Mat F,Mat A,const MatFactorInfo *info) 1318ab235cb6SHong Zhang { 1319a093e273SMatthew Knepley Mat B = F; 1320ab235cb6SHong Zhang Mat_Plapack *lu; 1321ab235cb6SHong Zhang PetscErrorCode ierr; 1322d0f46423SBarry Smith PetscInt M=A->rmap->N; 1323ab235cb6SHong Zhang MPI_Comm comm=((PetscObject)A)->comm,comm_2d; 1324ab235cb6SHong Zhang PetscMPIInt size; 1325ab235cb6SHong Zhang PetscInt ierror; 1326ab235cb6SHong Zhang 1327ab235cb6SHong Zhang PetscFunctionBegin; 1328ab235cb6SHong Zhang lu = (Mat_Plapack*)(B->spptr); 1329ab235cb6SHong Zhang 1330ab235cb6SHong Zhang /* Set default Plapack parameters */ 1331ab235cb6SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1332ab235cb6SHong Zhang lu->nprows = 1; lu->npcols = size; 1333ab235cb6SHong Zhang ierror = 0; 1334ab235cb6SHong Zhang lu->nb = M/size; 1335ab235cb6SHong Zhang if (M - lu->nb*size) lu->nb++; /* without cyclic distribution */ 1336ab235cb6SHong Zhang 1337ab235cb6SHong Zhang /* Set runtime options */ 1338ab235cb6SHong Zhang ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PLAPACK Options","Mat");CHKERRQ(ierr); 1339ab235cb6SHong Zhang ierr = PetscOptionsInt("-mat_plapack_nprows","row dimension of 2D processor mesh","None",lu->nprows,&lu->nprows,PETSC_NULL);CHKERRQ(ierr); 1340ab235cb6SHong Zhang ierr = PetscOptionsInt("-mat_plapack_npcols","column dimension of 2D processor mesh","None",lu->npcols,&lu->npcols,PETSC_NULL);CHKERRQ(ierr); 1341ab235cb6SHong Zhang 1342ab235cb6SHong Zhang ierr = PetscOptionsInt("-mat_plapack_nb","block size of template vector","None",lu->nb,&lu->nb,PETSC_NULL);CHKERRQ(ierr); 1343ab235cb6SHong Zhang ierr = PetscOptionsInt("-mat_plapack_ckerror","error checking flag","None",ierror,&ierror,PETSC_NULL);CHKERRQ(ierr); 1344ab235cb6SHong Zhang if (ierror){ 1345ab235cb6SHong Zhang PLA_Set_error_checking(ierror,PETSC_TRUE,PETSC_TRUE,PETSC_FALSE ); 1346ab235cb6SHong Zhang } else { 1347ab235cb6SHong Zhang PLA_Set_error_checking(ierror,PETSC_FALSE,PETSC_FALSE,PETSC_FALSE ); 1348ab235cb6SHong Zhang } 1349ab235cb6SHong Zhang lu->ierror = ierror; 1350ab235cb6SHong Zhang 1351ab235cb6SHong Zhang lu->nb_alg = 0; 1352ab235cb6SHong Zhang ierr = PetscOptionsInt("-mat_plapack_nb_alg","algorithmic block size","None",lu->nb_alg,&lu->nb_alg,PETSC_NULL);CHKERRQ(ierr); 1353ab235cb6SHong Zhang if (lu->nb_alg){ 1354ab235cb6SHong Zhang pla_Environ_set_nb_alg (PLA_OP_ALL_ALG,lu->nb_alg); 1355ab235cb6SHong Zhang } 1356ab235cb6SHong Zhang PetscOptionsEnd(); 1357ab235cb6SHong Zhang 1358ab235cb6SHong Zhang 1359ab235cb6SHong Zhang /* Create a 2D communicator */ 1360ab235cb6SHong Zhang PLA_Comm_1D_to_2D(comm,lu->nprows,lu->npcols,&comm_2d); 1361ab235cb6SHong Zhang lu->comm_2d = comm_2d; 1362ab235cb6SHong Zhang 1363ab235cb6SHong Zhang /* Create object distribution template */ 1364ab235cb6SHong Zhang lu->templ = NULL; 1365ab235cb6SHong Zhang PLA_Temp_create(lu->nb, 0, &lu->templ); 1366ab235cb6SHong Zhang 1367ab235cb6SHong Zhang /* Use suggested nb_alg if it is not provided by user */ 1368ab235cb6SHong Zhang if (lu->nb_alg == 0){ 1369ab235cb6SHong Zhang PLA_Environ_nb_alg(PLA_OP_PAN_PAN,lu->templ,&lu->nb_alg); 1370ab235cb6SHong Zhang pla_Environ_set_nb_alg(PLA_OP_ALL_ALG,lu->nb_alg); 1371ab235cb6SHong Zhang } 1372ab235cb6SHong Zhang 1373ab235cb6SHong Zhang /* Set the datatype */ 1374ab235cb6SHong Zhang #if defined(PETSC_USE_COMPLEX) 1375ab235cb6SHong Zhang lu->datatype = MPI_DOUBLE_COMPLEX; 1376ab235cb6SHong Zhang #else 1377ab235cb6SHong Zhang lu->datatype = MPI_DOUBLE; 1378ab235cb6SHong Zhang #endif 1379ab235cb6SHong Zhang 1380ab235cb6SHong Zhang lu->pla_solved = PETSC_FALSE; /* MatSolve_Plapack() is called yet */ 1381ab235cb6SHong Zhang lu->mstruct = DIFFERENT_NONZERO_PATTERN; 1382ab235cb6SHong Zhang lu->CleanUpPlapack = PETSC_TRUE; 1383ab235cb6SHong Zhang PetscFunctionReturn(0); 1384ab235cb6SHong Zhang } 1385ab235cb6SHong Zhang 1386b24902e0SBarry Smith /* Note the Petsc perm permutation is ignored */ 13872fbe02b9SBarry Smith #undef __FUNCT__ 13886c995c7dSHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_MPIDense" 13890481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_MPIDense(Mat F,Mat A,IS perm,const MatFactorInfo *info) 139001b82886SBarry Smith { 139101b82886SBarry Smith PetscErrorCode ierr; 1392b24902e0SBarry Smith PetscTruth issymmetric,set; 139301b82886SBarry Smith 139401b82886SBarry Smith PetscFunctionBegin; 1395b24902e0SBarry Smith ierr = MatIsSymmetricKnown(A,&set,&issymmetric); CHKERRQ(ierr); 1396b24902e0SBarry Smith if (!set || !issymmetric) SETERRQ(PETSC_ERR_USER,"Matrix must be set as MAT_SYMMETRIC for CholeskyFactor()"); 1397719d5645SBarry Smith ierr = MatFactorSymbolic_Plapack_Private(F,A,info);CHKERRQ(ierr); 1398719d5645SBarry Smith F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_MPIDense; 1399b24902e0SBarry Smith PetscFunctionReturn(0); 140001b82886SBarry Smith } 140101b82886SBarry Smith 1402b24902e0SBarry Smith /* Note the Petsc r and c permutations are ignored */ 1403b24902e0SBarry Smith #undef __FUNCT__ 14046c995c7dSHong Zhang #define __FUNCT__ "MatLUFactorSymbolic_MPIDense" 14050481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_MPIDense(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1406b24902e0SBarry Smith { 1407b24902e0SBarry Smith PetscErrorCode ierr; 1408d0f46423SBarry Smith PetscInt M = A->rmap->N; 1409b24902e0SBarry Smith Mat_Plapack *lu; 141001b82886SBarry Smith 1411b24902e0SBarry Smith PetscFunctionBegin; 1412a093e273SMatthew Knepley ierr = MatFactorSymbolic_Plapack_Private(F,A,info);CHKERRQ(ierr); 1413719d5645SBarry Smith lu = (Mat_Plapack*)F->spptr; 1414b24902e0SBarry Smith ierr = PLA_Mvector_create(MPI_INT,M,1,lu->templ,PLA_ALIGN_FIRST,&lu->pivots);CHKERRQ(ierr); 1415719d5645SBarry Smith F->ops->lufactornumeric = MatLUFactorNumeric_MPIDense; 141601b82886SBarry Smith PetscFunctionReturn(0); 141701b82886SBarry Smith } 141801b82886SBarry Smith 141901b82886SBarry Smith #undef __FUNCT__ 142085bd4cefSHong Zhang #define __FUNCT__ "MatGetFactor_mpidense_petsc" 142185bd4cefSHong Zhang PetscErrorCode MatGetFactor_mpidense_petsc(Mat A,MatFactorType ftype,Mat *F) 142201b82886SBarry Smith { 142301b82886SBarry Smith PetscErrorCode ierr; 14246c995c7dSHong Zhang Mat_Plapack *lu; 14256c995c7dSHong Zhang PetscMPIInt size; 142685bd4cefSHong Zhang PetscInt M=A->rmap->N; 142701b82886SBarry Smith 142801b82886SBarry Smith PetscFunctionBegin; 142901b82886SBarry Smith /* Create the factorization matrix */ 1430eae6fb2eSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,F);CHKERRQ(ierr); 1431d0f46423SBarry Smith ierr = MatSetSizes(*F,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1432eae6fb2eSBarry Smith ierr = MatSetType(*F,((PetscObject)A)->type_name);CHKERRQ(ierr); 143385bd4cefSHong Zhang ierr = PetscNewLog(*F,Mat_Plapack,&lu);CHKERRQ(ierr); 143485bd4cefSHong Zhang (*F)->spptr = (void*)lu; 143501b82886SBarry Smith 1436b24902e0SBarry Smith /* Set default Plapack parameters */ 14376c995c7dSHong Zhang ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1438b24902e0SBarry Smith lu->nb = M/size; 1439b24902e0SBarry Smith if (M - lu->nb*size) lu->nb++; /* without cyclic distribution */ 144001b82886SBarry Smith 1441b24902e0SBarry Smith /* Set runtime options */ 1442b24902e0SBarry Smith ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PLAPACK Options","Mat");CHKERRQ(ierr); 1443b24902e0SBarry Smith ierr = PetscOptionsInt("-mat_plapack_nb","block size of template vector","None",lu->nb,&lu->nb,PETSC_NULL);CHKERRQ(ierr); 1444b24902e0SBarry Smith PetscOptionsEnd(); 144501b82886SBarry Smith 1446b24902e0SBarry Smith /* Create object distribution template */ 1447b24902e0SBarry Smith lu->templ = NULL; 1448b24902e0SBarry Smith ierr = PLA_Temp_create(lu->nb, 0, &lu->templ);CHKERRQ(ierr); 1449b24902e0SBarry Smith 1450b24902e0SBarry Smith /* Set the datatype */ 1451b24902e0SBarry Smith #if defined(PETSC_USE_COMPLEX) 1452b24902e0SBarry Smith lu->datatype = MPI_DOUBLE_COMPLEX; 1453b24902e0SBarry Smith #else 1454b24902e0SBarry Smith lu->datatype = MPI_DOUBLE; 1455b24902e0SBarry Smith #endif 1456b24902e0SBarry Smith 1457d0f46423SBarry Smith ierr = PLA_Matrix_create(lu->datatype,M,A->cmap->N,lu->templ,PLA_ALIGN_FIRST,PLA_ALIGN_FIRST,&lu->A);CHKERRQ(ierr); 1458b24902e0SBarry Smith 1459b24902e0SBarry Smith 1460b24902e0SBarry Smith lu->pla_solved = PETSC_FALSE; /* MatSolve_Plapack() is called yet */ 1461b24902e0SBarry Smith 1462b24902e0SBarry Smith if (ftype == MAT_FACTOR_LU) { 1463b24902e0SBarry Smith (*F)->ops->lufactorsymbolic = MatLUFactorSymbolic_MPIDense; 1464b24902e0SBarry Smith } else if (ftype == MAT_FACTOR_CHOLESKY) { 14656c995c7dSHong Zhang (*F)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MPIDense; 1466b24902e0SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"No incomplete factorizations for dense matrices"); 1467719d5645SBarry Smith (*F)->factor = ftype; 146801b82886SBarry Smith PetscFunctionReturn(0); 146901b82886SBarry Smith } 147001b82886SBarry Smith #endif 147101b82886SBarry Smith 14728965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 147309dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 147409dc0095SBarry Smith MatGetRow_MPIDense, 147509dc0095SBarry Smith MatRestoreRow_MPIDense, 147609dc0095SBarry Smith MatMult_MPIDense, 147797304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 14787c922b88SBarry Smith MatMultTranspose_MPIDense, 14797c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 14808965ea79SLois Curfman McInnes 0, 148109dc0095SBarry Smith 0, 148209dc0095SBarry Smith 0, 148397304618SKris Buschelman /*10*/ 0, 148409dc0095SBarry Smith 0, 148509dc0095SBarry Smith 0, 148609dc0095SBarry Smith 0, 148709dc0095SBarry Smith MatTranspose_MPIDense, 148897304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 14896e4ee0c6SHong Zhang MatEqual_MPIDense, 149009dc0095SBarry Smith MatGetDiagonal_MPIDense, 14915b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 149209dc0095SBarry Smith MatNorm_MPIDense, 149397304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 149409dc0095SBarry Smith MatAssemblyEnd_MPIDense, 149509dc0095SBarry Smith 0, 149609dc0095SBarry Smith MatSetOption_MPIDense, 149709dc0095SBarry Smith MatZeroEntries_MPIDense, 149897304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 1499919b68f7SBarry Smith 0, 150001b82886SBarry Smith 0, 150101b82886SBarry Smith 0, 150201b82886SBarry Smith 0, 150397304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 1504273d9f13SBarry Smith 0, 150509dc0095SBarry Smith 0, 150609dc0095SBarry Smith MatGetArray_MPIDense, 150709dc0095SBarry Smith MatRestoreArray_MPIDense, 150897304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 150909dc0095SBarry Smith 0, 151009dc0095SBarry Smith 0, 151109dc0095SBarry Smith 0, 151209dc0095SBarry Smith 0, 151397304618SKris Buschelman /*40*/ 0, 15142ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 151509dc0095SBarry Smith 0, 151609dc0095SBarry Smith MatGetValues_MPIDense, 151709dc0095SBarry Smith 0, 151897304618SKris Buschelman /*45*/ 0, 151909dc0095SBarry Smith MatScale_MPIDense, 152009dc0095SBarry Smith 0, 152109dc0095SBarry Smith 0, 152209dc0095SBarry Smith 0, 1523521d7252SBarry Smith /*50*/ 0, 152409dc0095SBarry Smith 0, 152509dc0095SBarry Smith 0, 152609dc0095SBarry Smith 0, 152709dc0095SBarry Smith 0, 152897304618SKris Buschelman /*55*/ 0, 152909dc0095SBarry Smith 0, 153009dc0095SBarry Smith 0, 153109dc0095SBarry Smith 0, 153209dc0095SBarry Smith 0, 153397304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 1534b9b97703SBarry Smith MatDestroy_MPIDense, 1535b9b97703SBarry Smith MatView_MPIDense, 1536357abbc8SBarry Smith 0, 153797304618SKris Buschelman 0, 153897304618SKris Buschelman /*65*/ 0, 153997304618SKris Buschelman 0, 154097304618SKris Buschelman 0, 154197304618SKris Buschelman 0, 154297304618SKris Buschelman 0, 154397304618SKris Buschelman /*70*/ 0, 154497304618SKris Buschelman 0, 154597304618SKris Buschelman 0, 154697304618SKris Buschelman 0, 154797304618SKris Buschelman 0, 154897304618SKris Buschelman /*75*/ 0, 154997304618SKris Buschelman 0, 155097304618SKris Buschelman 0, 155197304618SKris Buschelman 0, 155297304618SKris Buschelman 0, 155397304618SKris Buschelman /*80*/ 0, 155497304618SKris Buschelman 0, 155597304618SKris Buschelman 0, 155697304618SKris Buschelman 0, 1557865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense, 1558865e5f61SKris Buschelman 0, 1559865e5f61SKris Buschelman 0, 1560865e5f61SKris Buschelman 0, 1561865e5f61SKris Buschelman 0, 1562865e5f61SKris Buschelman 0, 15632fbe02b9SBarry Smith /*90*/ 15642fbe02b9SBarry Smith #if defined(PETSC_HAVE_PLAPACK) 15652fbe02b9SBarry Smith MatMatMult_MPIDense_MPIDense, 15664ae313f4SHong Zhang MatMatMultSymbolic_MPIDense_MPIDense, 15672fbe02b9SBarry Smith MatMatMultNumeric_MPIDense_MPIDense, 15682fbe02b9SBarry Smith #else 1569865e5f61SKris Buschelman 0, 1570865e5f61SKris Buschelman 0, 1571865e5f61SKris Buschelman 0, 15722fbe02b9SBarry Smith #endif 15732fbe02b9SBarry Smith 0, 1574865e5f61SKris Buschelman /*95*/ 0, 1575865e5f61SKris Buschelman 0, 1576865e5f61SKris Buschelman 0, 1577865e5f61SKris Buschelman 0}; 15788965ea79SLois Curfman McInnes 1579273d9f13SBarry Smith EXTERN_C_BEGIN 15804a2ae208SSatish Balay #undef __FUNCT__ 1581a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 1582be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1583a23d5eceSKris Buschelman { 1584a23d5eceSKris Buschelman Mat_MPIDense *a; 1585dfbe8321SBarry Smith PetscErrorCode ierr; 1586a23d5eceSKris Buschelman 1587a23d5eceSKris Buschelman PetscFunctionBegin; 1588a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1589a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1590a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1591a23d5eceSKris Buschelman 1592a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 1593f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1594d0f46423SBarry Smith ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr); 15955c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 15965c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 159752e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 1598a23d5eceSKris Buschelman PetscFunctionReturn(0); 1599a23d5eceSKris Buschelman } 1600a23d5eceSKris Buschelman EXTERN_C_END 1601a23d5eceSKris Buschelman 16020bad9183SKris Buschelman /*MC 1603fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 16040bad9183SKris Buschelman 16050bad9183SKris Buschelman Options Database Keys: 16060bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 16070bad9183SKris Buschelman 16080bad9183SKris Buschelman Level: beginner 16090bad9183SKris Buschelman 16107878bbefSBarry Smith MATMPIDENSE matrices may use direct solvers (LU, Cholesky, and QR) 16117878bbefSBarry Smith for parallel dense matrices via the external package PLAPACK, if PLAPACK is installed 16127878bbefSBarry Smith (run config/configure.py with the option --download-plapack) 16137878bbefSBarry Smith 16147878bbefSBarry Smith 16157878bbefSBarry Smith Options Database Keys: 16167878bbefSBarry Smith . -mat_plapack_nprows <n> - number of rows in processor partition 16177878bbefSBarry Smith . -mat_plapack_npcols <n> - number of columns in processor partition 16187878bbefSBarry Smith . -mat_plapack_nb <n> - block size of template vector 16197878bbefSBarry Smith . -mat_plapack_nb_alg <n> - algorithmic block size 16207878bbefSBarry Smith - -mat_plapack_ckerror <n> - error checking flag 16217878bbefSBarry Smith 16227878bbefSBarry Smith .seealso: MatCreateMPIDense(), MATDENSE, MATSEQDENSE 16230bad9183SKris Buschelman M*/ 16240bad9183SKris Buschelman 1625a23d5eceSKris Buschelman EXTERN_C_BEGIN 1626a23d5eceSKris Buschelman #undef __FUNCT__ 16274a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1628be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat) 1629273d9f13SBarry Smith { 1630273d9f13SBarry Smith Mat_MPIDense *a; 1631dfbe8321SBarry Smith PetscErrorCode ierr; 1632273d9f13SBarry Smith 1633273d9f13SBarry Smith PetscFunctionBegin; 163438f2d2fdSLisandro Dalcin ierr = PetscNewLog(mat,Mat_MPIDense,&a);CHKERRQ(ierr); 1635b0a32e0cSBarry Smith mat->data = (void*)a; 1636273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1637273d9f13SBarry Smith mat->mapping = 0; 1638273d9f13SBarry Smith 1639273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 16407adad957SLisandro Dalcin ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&a->rank);CHKERRQ(ierr); 16417adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)mat)->comm,&a->size);CHKERRQ(ierr); 1642273d9f13SBarry Smith 1643d0f46423SBarry Smith mat->rmap->bs = mat->cmap->bs = 1; 1644d0f46423SBarry Smith ierr = PetscMapSetUp(mat->rmap);CHKERRQ(ierr); 1645d0f46423SBarry Smith ierr = PetscMapSetUp(mat->cmap);CHKERRQ(ierr); 1646d0f46423SBarry Smith a->nvec = mat->cmap->n; 1647273d9f13SBarry Smith 1648273d9f13SBarry Smith /* build cache for off array entries formed */ 1649273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 16507adad957SLisandro Dalcin ierr = MatStashCreate_Private(((PetscObject)mat)->comm,1,&mat->stash);CHKERRQ(ierr); 1651273d9f13SBarry Smith 1652273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1653273d9f13SBarry Smith a->lvec = 0; 1654273d9f13SBarry Smith a->Mvctx = 0; 1655273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1656273d9f13SBarry Smith 1657273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1658273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1659273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1660a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1661a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1662a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 16634ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C", 16644ae313f4SHong Zhang "MatMatMult_MPIAIJ_MPIDense", 16654ae313f4SHong Zhang MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr); 16664ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C", 16674ae313f4SHong Zhang "MatMatMultSymbolic_MPIAIJ_MPIDense", 16684ae313f4SHong Zhang MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr); 16694ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C", 16704ae313f4SHong Zhang "MatMatMultNumeric_MPIAIJ_MPIDense", 16714ae313f4SHong Zhang MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr); 16727878bbefSBarry Smith #if defined(PETSC_HAVE_PLAPACK) 167385bd4cefSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetFactor_mpidense_petsc_C", 167485bd4cefSHong Zhang "MatGetFactor_mpidense_petsc", 167585bd4cefSHong Zhang MatGetFactor_mpidense_petsc);CHKERRQ(ierr); 1676db4efbfdSBarry Smith #if defined(PETSC_HAVE_PLAPACK) 1677db4efbfdSBarry Smith ierr = PetscPLAPACKInitializePackage(((PetscObject)mat)->comm);CHKERRQ(ierr); 1678db4efbfdSBarry Smith #endif 1679db4efbfdSBarry Smith 16807878bbefSBarry Smith #endif 168138aed534SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr); 168201b82886SBarry Smith 1683273d9f13SBarry Smith PetscFunctionReturn(0); 1684273d9f13SBarry Smith } 1685273d9f13SBarry Smith EXTERN_C_END 1686273d9f13SBarry Smith 1687209238afSKris Buschelman /*MC 1688002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1689209238afSKris Buschelman 1690209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1691209238afSKris Buschelman and MATMPIDENSE otherwise. 1692209238afSKris Buschelman 1693209238afSKris Buschelman Options Database Keys: 1694209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1695209238afSKris Buschelman 1696209238afSKris Buschelman Level: beginner 1697209238afSKris Buschelman 169801b82886SBarry Smith 1699209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1700209238afSKris Buschelman M*/ 1701209238afSKris Buschelman 1702209238afSKris Buschelman EXTERN_C_BEGIN 1703209238afSKris Buschelman #undef __FUNCT__ 1704209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1705be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A) 1706dfbe8321SBarry Smith { 17076849ba73SBarry Smith PetscErrorCode ierr; 170813f74950SBarry Smith PetscMPIInt size; 1709209238afSKris Buschelman 1710209238afSKris Buschelman PetscFunctionBegin; 17117adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1712209238afSKris Buschelman if (size == 1) { 1713209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1714209238afSKris Buschelman } else { 1715209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1716209238afSKris Buschelman } 1717209238afSKris Buschelman PetscFunctionReturn(0); 1718209238afSKris Buschelman } 1719209238afSKris Buschelman EXTERN_C_END 1720209238afSKris Buschelman 17214a2ae208SSatish Balay #undef __FUNCT__ 17224a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1723273d9f13SBarry Smith /*@C 1724273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1725273d9f13SBarry Smith 1726273d9f13SBarry Smith Not collective 1727273d9f13SBarry Smith 1728273d9f13SBarry Smith Input Parameters: 1729273d9f13SBarry Smith . A - the matrix 1730273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1731273d9f13SBarry Smith to control all matrix memory allocation. 1732273d9f13SBarry Smith 1733273d9f13SBarry Smith Notes: 1734273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1735273d9f13SBarry Smith storage by columns. 1736273d9f13SBarry Smith 1737273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1738273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1739273d9f13SBarry Smith set data=PETSC_NULL. 1740273d9f13SBarry Smith 1741273d9f13SBarry Smith Level: intermediate 1742273d9f13SBarry Smith 1743273d9f13SBarry Smith .keywords: matrix,dense, parallel 1744273d9f13SBarry Smith 1745273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1746273d9f13SBarry Smith @*/ 1747be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1748273d9f13SBarry Smith { 1749dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar *); 1750273d9f13SBarry Smith 1751273d9f13SBarry Smith PetscFunctionBegin; 1752565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1753a23d5eceSKris Buschelman if (f) { 1754a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1755a23d5eceSKris Buschelman } 1756273d9f13SBarry Smith PetscFunctionReturn(0); 1757273d9f13SBarry Smith } 1758273d9f13SBarry Smith 17594a2ae208SSatish Balay #undef __FUNCT__ 17604a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 17618965ea79SLois Curfman McInnes /*@C 176239ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 17638965ea79SLois Curfman McInnes 1764db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1765db81eaa0SLois Curfman McInnes 17668965ea79SLois Curfman McInnes Input Parameters: 1767db81eaa0SLois Curfman McInnes + comm - MPI communicator 17688965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1769db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 17708965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1771db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 17727f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1773dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 17748965ea79SLois Curfman McInnes 17758965ea79SLois Curfman McInnes Output Parameter: 1776477f1c0bSLois Curfman McInnes . A - the matrix 17778965ea79SLois Curfman McInnes 1778b259b22eSLois Curfman McInnes Notes: 177939ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 178039ddd567SLois Curfman McInnes storage by columns. 17818965ea79SLois Curfman McInnes 178218f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 178318f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 17847f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 178518f449edSLois Curfman McInnes 17868965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 17878965ea79SLois Curfman McInnes (possibly both). 17888965ea79SLois Curfman McInnes 1789027ccd11SLois Curfman McInnes Level: intermediate 1790027ccd11SLois Curfman McInnes 179139ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 17928965ea79SLois Curfman McInnes 179339ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 17948965ea79SLois Curfman McInnes @*/ 1795be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 17968965ea79SLois Curfman McInnes { 17976849ba73SBarry Smith PetscErrorCode ierr; 179813f74950SBarry Smith PetscMPIInt size; 17998965ea79SLois Curfman McInnes 18003a40ed3dSBarry Smith PetscFunctionBegin; 1801f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1802f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1803273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1804273d9f13SBarry Smith if (size > 1) { 1805273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1806273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1807273d9f13SBarry Smith } else { 1808273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1809273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 18108c469469SLois Curfman McInnes } 18113a40ed3dSBarry Smith PetscFunctionReturn(0); 18128965ea79SLois Curfman McInnes } 18138965ea79SLois Curfman McInnes 18144a2ae208SSatish Balay #undef __FUNCT__ 18154a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 18166849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 18178965ea79SLois Curfman McInnes { 18188965ea79SLois Curfman McInnes Mat mat; 18193501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1820dfbe8321SBarry Smith PetscErrorCode ierr; 18218965ea79SLois Curfman McInnes 18223a40ed3dSBarry Smith PetscFunctionBegin; 18238965ea79SLois Curfman McInnes *newmat = 0; 18247adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&mat);CHKERRQ(ierr); 1825d0f46423SBarry Smith ierr = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 18267adad957SLisandro Dalcin ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr); 1827834f8fabSBarry Smith a = (Mat_MPIDense*)mat->data; 1828e04c1aa4SHong Zhang ierr = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 18293501a2bdSLois Curfman McInnes mat->factor = A->factor; 1830c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1831273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 18328965ea79SLois Curfman McInnes 1833d0f46423SBarry Smith mat->rmap->rstart = A->rmap->rstart; 1834d0f46423SBarry Smith mat->rmap->rend = A->rmap->rend; 18358965ea79SLois Curfman McInnes a->size = oldmat->size; 18368965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1837e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1838b9b97703SBarry Smith a->nvec = oldmat->nvec; 18393782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1840e04c1aa4SHong Zhang 1841d0f46423SBarry Smith ierr = PetscMemcpy(mat->rmap->range,A->rmap->range,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr); 1842d0f46423SBarry Smith ierr = PetscMemcpy(mat->cmap->range,A->cmap->range,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr); 18437adad957SLisandro Dalcin ierr = MatStashCreate_Private(((PetscObject)A)->comm,1,&mat->stash);CHKERRQ(ierr); 18448965ea79SLois Curfman McInnes 1845329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 18465609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 184752e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 184801b82886SBarry Smith 18498965ea79SLois Curfman McInnes *newmat = mat; 18503a40ed3dSBarry Smith PetscFunctionReturn(0); 18518965ea79SLois Curfman McInnes } 18528965ea79SLois Curfman McInnes 1853e090d566SSatish Balay #include "petscsys.h" 18548965ea79SLois Curfman McInnes 18554a2ae208SSatish Balay #undef __FUNCT__ 18564a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 1857a313700dSBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, const MatType type,Mat *newmat) 185890ace30eSBarry Smith { 18596849ba73SBarry Smith PetscErrorCode ierr; 186032dcc486SBarry Smith PetscMPIInt rank,size; 186113f74950SBarry Smith PetscInt *rowners,i,m,nz,j; 186287828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 186390ace30eSBarry Smith MPI_Status status; 186490ace30eSBarry Smith 18653a40ed3dSBarry Smith PetscFunctionBegin; 1866d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1867d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 186890ace30eSBarry Smith 186990ace30eSBarry Smith /* determine ownership of all rows */ 187090ace30eSBarry Smith m = M/size + ((M % size) > rank); 187113f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 187213f74950SBarry Smith ierr = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 187390ace30eSBarry Smith rowners[0] = 0; 187490ace30eSBarry Smith for (i=2; i<=size; i++) { 187590ace30eSBarry Smith rowners[i] += rowners[i-1]; 187690ace30eSBarry Smith } 187790ace30eSBarry Smith 1878f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1879f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1880be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1881878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 188290ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 188390ace30eSBarry Smith 188490ace30eSBarry Smith if (!rank) { 188587828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 188690ace30eSBarry Smith 188790ace30eSBarry Smith /* read in my part of the matrix numerical values */ 18880752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 188990ace30eSBarry Smith 189090ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 189190ace30eSBarry Smith vals_ptr = vals; 189290ace30eSBarry Smith for (i=0; i<m; i++) { 189390ace30eSBarry Smith for (j=0; j<N; j++) { 189490ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 189590ace30eSBarry Smith } 189690ace30eSBarry Smith } 189790ace30eSBarry Smith 189890ace30eSBarry Smith /* read in other processors and ship out */ 189990ace30eSBarry Smith for (i=1; i<size; i++) { 190090ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 19010752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 19027adad957SLisandro Dalcin ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(*newmat))->tag,comm);CHKERRQ(ierr); 190390ace30eSBarry Smith } 19043a40ed3dSBarry Smith } else { 190590ace30eSBarry Smith /* receive numeric values */ 190687828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 190790ace30eSBarry Smith 190890ace30eSBarry Smith /* receive message of values*/ 19097adad957SLisandro Dalcin ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(*newmat))->tag,comm,&status);CHKERRQ(ierr); 191090ace30eSBarry Smith 191190ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 191290ace30eSBarry Smith vals_ptr = vals; 191390ace30eSBarry Smith for (i=0; i<m; i++) { 191490ace30eSBarry Smith for (j=0; j<N; j++) { 191590ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 191690ace30eSBarry Smith } 191790ace30eSBarry Smith } 191890ace30eSBarry Smith } 1919606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1920606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 19216d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19226d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19233a40ed3dSBarry Smith PetscFunctionReturn(0); 192490ace30eSBarry Smith } 192590ace30eSBarry Smith 19264a2ae208SSatish Balay #undef __FUNCT__ 19274a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1928a313700dSBarry Smith PetscErrorCode MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat) 19298965ea79SLois Curfman McInnes { 19308965ea79SLois Curfman McInnes Mat A; 193187828ca2SBarry Smith PetscScalar *vals,*svals; 193219bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 19338965ea79SLois Curfman McInnes MPI_Status status; 193413f74950SBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz; 193513f74950SBarry Smith PetscInt header[4],*rowlengths = 0,M,N,*cols; 193613f74950SBarry Smith PetscInt *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols; 193713f74950SBarry Smith PetscInt i,nz,j,rstart,rend; 193813f74950SBarry Smith int fd; 19396849ba73SBarry Smith PetscErrorCode ierr; 19408965ea79SLois Curfman McInnes 19413a40ed3dSBarry Smith PetscFunctionBegin; 1942d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1943d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 19448965ea79SLois Curfman McInnes if (!rank) { 1945b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 19460752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1947552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 19488965ea79SLois Curfman McInnes } 19498965ea79SLois Curfman McInnes 195013f74950SBarry Smith ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 195190ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 195290ace30eSBarry Smith 195390ace30eSBarry Smith /* 195490ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 195590ace30eSBarry Smith */ 195690ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1957be5d1d56SKris Buschelman ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 19583a40ed3dSBarry Smith PetscFunctionReturn(0); 195990ace30eSBarry Smith } 196090ace30eSBarry Smith 19618965ea79SLois Curfman McInnes /* determine ownership of all rows */ 1962e44c0bd4SBarry Smith m = PetscMPIIntCast(M/size + ((M % size) > rank)); 1963e44c0bd4SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscMPIInt),&rowners);CHKERRQ(ierr); 1964ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 19658965ea79SLois Curfman McInnes rowners[0] = 0; 19668965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 19678965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 19688965ea79SLois Curfman McInnes } 19698965ea79SLois Curfman McInnes rstart = rowners[rank]; 19708965ea79SLois Curfman McInnes rend = rowners[rank+1]; 19718965ea79SLois Curfman McInnes 19728965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 197313f74950SBarry Smith ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr); 19748965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 19758965ea79SLois Curfman McInnes if (!rank) { 197613f74950SBarry Smith ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 19770752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 197813f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr); 19798965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 19801466f1e1SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 1981606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1982ca161407SBarry Smith } else { 19831466f1e1SBarry Smith ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 19848965ea79SLois Curfman McInnes } 19858965ea79SLois Curfman McInnes 19868965ea79SLois Curfman McInnes if (!rank) { 19878965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 198813f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr); 198913f74950SBarry Smith ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 19908965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 19918965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 19928965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 19938965ea79SLois Curfman McInnes } 19948965ea79SLois Curfman McInnes } 1995606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 19968965ea79SLois Curfman McInnes 19978965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 19988965ea79SLois Curfman McInnes maxnz = 0; 19998965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 20000452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 20018965ea79SLois Curfman McInnes } 200213f74950SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr); 20038965ea79SLois Curfman McInnes 20048965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 20058965ea79SLois Curfman McInnes nz = procsnz[0]; 200613f74950SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 20070752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 20088965ea79SLois Curfman McInnes 20098965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 20108965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 20118965ea79SLois Curfman McInnes nz = procsnz[i]; 20120752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 201313f74950SBarry Smith ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 20148965ea79SLois Curfman McInnes } 2015606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 20163a40ed3dSBarry Smith } else { 20178965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 20188965ea79SLois Curfman McInnes nz = 0; 20198965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 20208965ea79SLois Curfman McInnes nz += ourlens[i]; 20218965ea79SLois Curfman McInnes } 202213f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 20238965ea79SLois Curfman McInnes 20248965ea79SLois Curfman McInnes /* receive message of column indices*/ 202513f74950SBarry Smith ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 202613f74950SBarry Smith ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 202729bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 20288965ea79SLois Curfman McInnes } 20298965ea79SLois Curfman McInnes 20308965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 203113f74950SBarry Smith ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr); 20328965ea79SLois Curfman McInnes jj = 0; 20338965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 20348965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 20358965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 20368965ea79SLois Curfman McInnes jj++; 20378965ea79SLois Curfman McInnes } 20388965ea79SLois Curfman McInnes } 20398965ea79SLois Curfman McInnes 20408965ea79SLois Curfman McInnes /* create our matrix */ 20418965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 20428965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 20438965ea79SLois Curfman McInnes } 2044f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 2045f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 2046be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 2047878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 20488965ea79SLois Curfman McInnes A = *newmat; 20498965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 20508965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 20518965ea79SLois Curfman McInnes } 20528965ea79SLois Curfman McInnes 20538965ea79SLois Curfman McInnes if (!rank) { 205487828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 20558965ea79SLois Curfman McInnes 20568965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 20578965ea79SLois Curfman McInnes nz = procsnz[0]; 20580752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 20598965ea79SLois Curfman McInnes 20608965ea79SLois Curfman McInnes /* insert into matrix */ 20618965ea79SLois Curfman McInnes jj = rstart; 20628965ea79SLois Curfman McInnes smycols = mycols; 20638965ea79SLois Curfman McInnes svals = vals; 20648965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 20658965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 20668965ea79SLois Curfman McInnes smycols += ourlens[i]; 20678965ea79SLois Curfman McInnes svals += ourlens[i]; 20688965ea79SLois Curfman McInnes jj++; 20698965ea79SLois Curfman McInnes } 20708965ea79SLois Curfman McInnes 20718965ea79SLois Curfman McInnes /* read in other processors and ship out */ 20728965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 20738965ea79SLois Curfman McInnes nz = procsnz[i]; 20740752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 20757adad957SLisandro Dalcin ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)A)->tag,comm);CHKERRQ(ierr); 20768965ea79SLois Curfman McInnes } 2077606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 20783a40ed3dSBarry Smith } else { 20798965ea79SLois Curfman McInnes /* receive numeric values */ 208084d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 20818965ea79SLois Curfman McInnes 20828965ea79SLois Curfman McInnes /* receive message of values*/ 20837adad957SLisandro Dalcin ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)A)->tag,comm,&status);CHKERRQ(ierr); 2084ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 208529bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 20868965ea79SLois Curfman McInnes 20878965ea79SLois Curfman McInnes /* insert into matrix */ 20888965ea79SLois Curfman McInnes jj = rstart; 20898965ea79SLois Curfman McInnes smycols = mycols; 20908965ea79SLois Curfman McInnes svals = vals; 20918965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 20928965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 20938965ea79SLois Curfman McInnes smycols += ourlens[i]; 20948965ea79SLois Curfman McInnes svals += ourlens[i]; 20958965ea79SLois Curfman McInnes jj++; 20968965ea79SLois Curfman McInnes } 20978965ea79SLois Curfman McInnes } 2098606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 2099606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 2100606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 2101606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 21028965ea79SLois Curfman McInnes 21036d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 21046d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 21053a40ed3dSBarry Smith PetscFunctionReturn(0); 21068965ea79SLois Curfman McInnes } 210790ace30eSBarry Smith 21086e4ee0c6SHong Zhang #undef __FUNCT__ 21096e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense" 21106e4ee0c6SHong Zhang PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscTruth *flag) 21116e4ee0c6SHong Zhang { 21126e4ee0c6SHong Zhang Mat_MPIDense *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data; 21136e4ee0c6SHong Zhang Mat a,b; 21146e4ee0c6SHong Zhang PetscTruth flg; 21156e4ee0c6SHong Zhang PetscErrorCode ierr; 211690ace30eSBarry Smith 21176e4ee0c6SHong Zhang PetscFunctionBegin; 21186e4ee0c6SHong Zhang a = matA->A; 21196e4ee0c6SHong Zhang b = matB->A; 21206e4ee0c6SHong Zhang ierr = MatEqual(a,b,&flg);CHKERRQ(ierr); 21217adad957SLisandro Dalcin ierr = MPI_Allreduce(&flg,flag,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr); 21226e4ee0c6SHong Zhang PetscFunctionReturn(0); 21236e4ee0c6SHong Zhang } 212490ace30eSBarry Smith 212509d27a7eSBarry Smith #if defined(PETSC_HAVE_PLAPACK) 212609d27a7eSBarry Smith 212709d27a7eSBarry Smith #undef __FUNCT__ 212809d27a7eSBarry Smith #define __FUNCT__ "PetscPLAPACKFinalizePackage" 212909d27a7eSBarry Smith /*@C 213009d27a7eSBarry Smith PetscPLAPACKFinalizePackage - This function destroys everything in the Petsc interface to PLAPACK. 213109d27a7eSBarry Smith Level: developer 213209d27a7eSBarry Smith 213309d27a7eSBarry Smith .keywords: Petsc, destroy, package, PLAPACK 213409d27a7eSBarry Smith .seealso: PetscFinalize() 213509d27a7eSBarry Smith @*/ 213609d27a7eSBarry Smith PetscErrorCode PETSC_DLLEXPORT PetscPLAPACKFinalizePackage(void) 213709d27a7eSBarry Smith { 213809d27a7eSBarry Smith PetscErrorCode ierr; 213909d27a7eSBarry Smith 214009d27a7eSBarry Smith PetscFunctionBegin; 214109d27a7eSBarry Smith ierr = PLA_Finalize();CHKERRQ(ierr); 214209d27a7eSBarry Smith PetscFunctionReturn(0); 214309d27a7eSBarry Smith } 214409d27a7eSBarry Smith 214509d27a7eSBarry Smith #undef __FUNCT__ 214609d27a7eSBarry Smith #define __FUNCT__ "PetscPLAPACKInitializePackage" 214709d27a7eSBarry Smith /*@C 214809d27a7eSBarry Smith PetscPLAPACKInitializePackage - This function initializes everything in the Petsc interface to PLAPACK. It is 2149db4efbfdSBarry Smith called from MatCreate_MPIDense() the first time an MPI dense matrix is called. 215009d27a7eSBarry Smith 215109d27a7eSBarry Smith Input Parameter: 2152db4efbfdSBarry Smith . comm - the communicator the matrix lives on 215309d27a7eSBarry Smith 215409d27a7eSBarry Smith Level: developer 215509d27a7eSBarry Smith 2156db4efbfdSBarry Smith Notes: PLAPACK does not have a good fit with MPI communicators; all (parallel) PLAPACK objects have to live in the 2157db4efbfdSBarry Smith same communicator (because there is some global state that is initialized and used for all matrices). In addition if 2158db4efbfdSBarry Smith PLAPACK is initialized (that is the initial matrices created) are on subcommunicators of MPI_COMM_WORLD, these subcommunicators 2159db4efbfdSBarry Smith cannot overlap. 2160db4efbfdSBarry Smith 216109d27a7eSBarry Smith .keywords: Petsc, initialize, package, PLAPACK 216209d27a7eSBarry Smith .seealso: PetscInitializePackage(), PetscInitialize() 216309d27a7eSBarry Smith @*/ 2164db4efbfdSBarry Smith PetscErrorCode PETSC_DLLEXPORT PetscPLAPACKInitializePackage(MPI_Comm comm) 216509d27a7eSBarry Smith { 2166eae6fb2eSBarry Smith PetscMPIInt size; 216709d27a7eSBarry Smith PetscErrorCode ierr; 216809d27a7eSBarry Smith 216909d27a7eSBarry Smith PetscFunctionBegin; 217009d27a7eSBarry Smith if (!PLA_Initialized(PETSC_NULL)) { 217109d27a7eSBarry Smith 217209d27a7eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 217304fea9ffSBarry Smith Plapack_nprows = 1; 217404fea9ffSBarry Smith Plapack_npcols = size; 217509d27a7eSBarry Smith 2176aeccfd6fSBarry Smith ierr = PetscOptionsBegin(comm,PETSC_NULL,"PLAPACK Options","Mat");CHKERRQ(ierr); 2177aeccfd6fSBarry Smith ierr = PetscOptionsInt("-plapack_nprows","row dimension of 2D processor mesh","None",Plapack_nprows,&Plapack_nprows,PETSC_NULL);CHKERRQ(ierr); 2178aeccfd6fSBarry Smith ierr = PetscOptionsInt("-plapack_npcols","column dimension of 2D processor mesh","None",Plapack_npcols,&Plapack_npcols,PETSC_NULL);CHKERRQ(ierr); 217979b0a62dSBarry Smith #if defined(PETSC_USE_DEBUG) 218079b0a62dSBarry Smith Plapack_ierror = 3; 218179b0a62dSBarry Smith #else 218279b0a62dSBarry Smith Plapack_ierror = 0; 218379b0a62dSBarry Smith #endif 2184eae6fb2eSBarry Smith ierr = PetscOptionsInt("-plapack_ckerror","error checking flag","None",Plapack_ierror,&Plapack_ierror,PETSC_NULL);CHKERRQ(ierr); 2185eae6fb2eSBarry Smith if (Plapack_ierror){ 2186eae6fb2eSBarry Smith ierr = PLA_Set_error_checking(Plapack_ierror,PETSC_TRUE,PETSC_TRUE,PETSC_FALSE );CHKERRQ(ierr); 2187aeccfd6fSBarry Smith } else { 2188eae6fb2eSBarry Smith ierr = PLA_Set_error_checking(Plapack_ierror,PETSC_FALSE,PETSC_FALSE,PETSC_FALSE );CHKERRQ(ierr); 2189aeccfd6fSBarry Smith } 2190aeccfd6fSBarry Smith 2191eae6fb2eSBarry Smith Plapack_nb_alg = 0; 219279b0a62dSBarry Smith ierr = PetscOptionsInt("-plapack_nb_alg","algorithmic block size","None",Plapack_nb_alg,&Plapack_nb_alg,PETSC_NULL);CHKERRQ(ierr); 2193eae6fb2eSBarry Smith if (Plapack_nb_alg) { 2194eae6fb2eSBarry Smith ierr = pla_Environ_set_nb_alg (PLA_OP_ALL_ALG,Plapack_nb_alg);CHKERRQ(ierr); 2195aeccfd6fSBarry Smith } 2196aeccfd6fSBarry Smith PetscOptionsEnd(); 2197aeccfd6fSBarry Smith 219804fea9ffSBarry Smith ierr = PLA_Comm_1D_to_2D(comm,Plapack_nprows,Plapack_npcols,&Plapack_comm_2d);CHKERRQ(ierr); 219904fea9ffSBarry Smith ierr = PLA_Init(Plapack_comm_2d);CHKERRQ(ierr); 220009d27a7eSBarry Smith ierr = PetscRegisterFinalize(PetscPLAPACKFinalizePackage);CHKERRQ(ierr); 220109d27a7eSBarry Smith } 220209d27a7eSBarry Smith PetscFunctionReturn(0); 220309d27a7eSBarry Smith } 220490ace30eSBarry Smith 220590ace30eSBarry Smith 220609d27a7eSBarry Smith #endif 2207