1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3ed3cc1f0SBarry Smith /* 4ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 5ed3cc1f0SBarry Smith */ 6ed3cc1f0SBarry Smith 789ac3759SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" /*I "petscmat.h" I*/ 88965ea79SLois Curfman McInnes 9ba8c8a56SBarry Smith #undef __FUNCT__ 10ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix" 11ab92ecdeSBarry Smith /*@ 12ab92ecdeSBarry Smith 13ab92ecdeSBarry Smith MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential 14ab92ecdeSBarry Smith matrix that represents the operator. For sequential matrices it returns itself. 15ab92ecdeSBarry Smith 16ab92ecdeSBarry Smith Input Parameter: 17ab92ecdeSBarry Smith . A - the Seq or MPI dense matrix 18ab92ecdeSBarry Smith 19ab92ecdeSBarry Smith Output Parameter: 20ab92ecdeSBarry Smith . B - the inner matrix 21ab92ecdeSBarry Smith 22*8e6c10adSSatish Balay Level: intermediate 23*8e6c10adSSatish Balay 24ab92ecdeSBarry Smith @*/ 25ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B) 26ab92ecdeSBarry Smith { 27ab92ecdeSBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 28ab92ecdeSBarry Smith PetscErrorCode ierr; 29ab92ecdeSBarry Smith PetscTruth flg; 30ab92ecdeSBarry Smith PetscMPIInt size; 31ab92ecdeSBarry Smith 32ab92ecdeSBarry Smith PetscFunctionBegin; 33ab92ecdeSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr); 34ab92ecdeSBarry Smith if (!flg) { /* this check sucks! */ 35ab92ecdeSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATDENSE,&flg);CHKERRQ(ierr); 36ab92ecdeSBarry Smith if (flg) { 37ab92ecdeSBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 38ab92ecdeSBarry Smith if (size == 1) flg = PETSC_FALSE; 39ab92ecdeSBarry Smith } 40ab92ecdeSBarry Smith } 41ab92ecdeSBarry Smith if (flg) { 42ab92ecdeSBarry Smith *B = mat->A; 43ab92ecdeSBarry Smith } else { 44ab92ecdeSBarry Smith *B = A; 45ab92ecdeSBarry Smith } 46ab92ecdeSBarry Smith PetscFunctionReturn(0); 47ab92ecdeSBarry Smith } 48ab92ecdeSBarry Smith 49ab92ecdeSBarry Smith #undef __FUNCT__ 50ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense" 51ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 52ba8c8a56SBarry Smith { 53ba8c8a56SBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 54ba8c8a56SBarry Smith PetscErrorCode ierr; 55ba8c8a56SBarry Smith PetscInt lrow,rstart = mat->rstart,rend = mat->rend; 56ba8c8a56SBarry Smith 57ba8c8a56SBarry Smith PetscFunctionBegin; 58ba8c8a56SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 59ba8c8a56SBarry Smith lrow = row - rstart; 60ba8c8a56SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr); 61ba8c8a56SBarry Smith PetscFunctionReturn(0); 62ba8c8a56SBarry Smith } 63ba8c8a56SBarry Smith 64ba8c8a56SBarry Smith #undef __FUNCT__ 65ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense" 66ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 67ba8c8a56SBarry Smith { 68ba8c8a56SBarry Smith PetscErrorCode ierr; 69ba8c8a56SBarry Smith 70ba8c8a56SBarry Smith PetscFunctionBegin; 71ba8c8a56SBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 72ba8c8a56SBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 73ba8c8a56SBarry Smith PetscFunctionReturn(0); 74ba8c8a56SBarry Smith } 75ba8c8a56SBarry Smith 760de54da6SSatish Balay EXTERN_C_BEGIN 774a2ae208SSatish Balay #undef __FUNCT__ 784a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 79be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 800de54da6SSatish Balay { 810de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 826849ba73SBarry Smith PetscErrorCode ierr; 8313f74950SBarry Smith PetscInt m = A->m,rstart = mdn->rstart; 8487828ca2SBarry Smith PetscScalar *array; 850de54da6SSatish Balay MPI_Comm comm; 860de54da6SSatish Balay 870de54da6SSatish Balay PetscFunctionBegin; 88273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 890de54da6SSatish Balay 900de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 910de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 920de54da6SSatish Balay 930de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 940de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 95f69a0ea3SMatthew Knepley ierr = MatCreate(comm,B);CHKERRQ(ierr); 96f69a0ea3SMatthew Knepley ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr); 975c5985e7SKris Buschelman ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr); 985c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr); 990de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 1000de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1010de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1020de54da6SSatish Balay 1030de54da6SSatish Balay *iscopy = PETSC_TRUE; 1040de54da6SSatish Balay PetscFunctionReturn(0); 1050de54da6SSatish Balay } 1060de54da6SSatish Balay EXTERN_C_END 1070de54da6SSatish Balay 1084a2ae208SSatish Balay #undef __FUNCT__ 1094a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 11013f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) 1118965ea79SLois Curfman McInnes { 11239b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 113dfbe8321SBarry Smith PetscErrorCode ierr; 11413f74950SBarry Smith PetscInt i,j,rstart = A->rstart,rend = A->rend,row; 115273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 1168965ea79SLois Curfman McInnes 1173a40ed3dSBarry Smith PetscFunctionBegin; 1188965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 1195ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 120273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 1218965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 1228965ea79SLois Curfman McInnes row = idxm[i] - rstart; 12339b7565bSBarry Smith if (roworiented) { 12439b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 1253a40ed3dSBarry Smith } else { 1268965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 1275ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 128273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 12939b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 13039b7565bSBarry Smith } 1318965ea79SLois Curfman McInnes } 1323a40ed3dSBarry Smith } else { 1333782ba37SSatish Balay if (!A->donotstash) { 13439b7565bSBarry Smith if (roworiented) { 1358798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 136d36fbae8SSatish Balay } else { 1378798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 13839b7565bSBarry Smith } 139b49de8d1SLois Curfman McInnes } 140b49de8d1SLois Curfman McInnes } 1413782ba37SSatish Balay } 1423a40ed3dSBarry Smith PetscFunctionReturn(0); 143b49de8d1SLois Curfman McInnes } 144b49de8d1SLois Curfman McInnes 1454a2ae208SSatish Balay #undef __FUNCT__ 1464a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 14713f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 148b49de8d1SLois Curfman McInnes { 149b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 150dfbe8321SBarry Smith PetscErrorCode ierr; 15113f74950SBarry Smith PetscInt i,j,rstart = mdn->rstart,rend = mdn->rend,row; 152b49de8d1SLois Curfman McInnes 1533a40ed3dSBarry Smith PetscFunctionBegin; 154b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 15529bbc08cSBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 156273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 157b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 158b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 159b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 16029bbc08cSBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 161273d9f13SBarry Smith if (idxn[j] >= mat->N) { 16229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 163a8c6a408SBarry Smith } 164b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 165b49de8d1SLois Curfman McInnes } 166a8c6a408SBarry Smith } else { 16729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 1688965ea79SLois Curfman McInnes } 1698965ea79SLois Curfman McInnes } 1703a40ed3dSBarry Smith PetscFunctionReturn(0); 1718965ea79SLois Curfman McInnes } 1728965ea79SLois Curfman McInnes 1734a2ae208SSatish Balay #undef __FUNCT__ 1744a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 175dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 176ff14e315SSatish Balay { 177ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 178dfbe8321SBarry Smith PetscErrorCode ierr; 179ff14e315SSatish Balay 1803a40ed3dSBarry Smith PetscFunctionBegin; 181ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1823a40ed3dSBarry Smith PetscFunctionReturn(0); 183ff14e315SSatish Balay } 184ff14e315SSatish Balay 1854a2ae208SSatish Balay #undef __FUNCT__ 1864a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 18713f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 188ca3fa75bSLois Curfman McInnes { 189ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 190ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 1916849ba73SBarry Smith PetscErrorCode ierr; 19213f74950SBarry Smith PetscInt i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 19387828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 194ca3fa75bSLois Curfman McInnes Mat newmat; 195ca3fa75bSLois Curfman McInnes 196ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 197ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 198ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 199b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 200b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 201ca3fa75bSLois Curfman McInnes 202ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 2037eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 2047eba5e9cSLois Curfman McInnes original matrix! */ 205ca3fa75bSLois Curfman McInnes 206ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 2077eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 208ca3fa75bSLois Curfman McInnes 209ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 210ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 21129bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 2127eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 213ca3fa75bSLois Curfman McInnes newmat = *B; 214ca3fa75bSLois Curfman McInnes } else { 215ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 216f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr); 217f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr); 218878740d9SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 219878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 220ca3fa75bSLois Curfman McInnes } 221ca3fa75bSLois Curfman McInnes 222ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 223ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 224ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 225ca3fa75bSLois Curfman McInnes 226ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 227ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 228ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 2297eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 230ca3fa75bSLois Curfman McInnes } 231ca3fa75bSLois Curfman McInnes } 232ca3fa75bSLois Curfman McInnes 233ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 234ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 235ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 236ca3fa75bSLois Curfman McInnes 237ca3fa75bSLois Curfman McInnes /* Free work space */ 238ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 239ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 240ca3fa75bSLois Curfman McInnes *B = newmat; 241ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 242ca3fa75bSLois Curfman McInnes } 243ca3fa75bSLois Curfman McInnes 2444a2ae208SSatish Balay #undef __FUNCT__ 2454a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 246dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 247ff14e315SSatish Balay { 2483a40ed3dSBarry Smith PetscFunctionBegin; 2493a40ed3dSBarry Smith PetscFunctionReturn(0); 250ff14e315SSatish Balay } 251ff14e315SSatish Balay 2524a2ae208SSatish Balay #undef __FUNCT__ 2534a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 254dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2558965ea79SLois Curfman McInnes { 25639ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 2578965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 258dfbe8321SBarry Smith PetscErrorCode ierr; 25913f74950SBarry Smith PetscInt nstash,reallocs; 2608965ea79SLois Curfman McInnes InsertMode addv; 2618965ea79SLois Curfman McInnes 2623a40ed3dSBarry Smith PetscFunctionBegin; 2638965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 264ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 2657056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 26629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 2678965ea79SLois Curfman McInnes } 268e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2698965ea79SLois Curfman McInnes 2708798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 2718798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 27263ba0a88SBarry Smith ierr = PetscLogInfo((mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs));CHKERRQ(ierr); 2733a40ed3dSBarry Smith PetscFunctionReturn(0); 2748965ea79SLois Curfman McInnes } 2758965ea79SLois Curfman McInnes 2764a2ae208SSatish Balay #undef __FUNCT__ 2774a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 278dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2798965ea79SLois Curfman McInnes { 28039ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2816849ba73SBarry Smith PetscErrorCode ierr; 28213f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 28313f74950SBarry Smith PetscMPIInt n; 28487828ca2SBarry Smith PetscScalar *val; 285e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2868965ea79SLois Curfman McInnes 2873a40ed3dSBarry Smith PetscFunctionBegin; 2888965ea79SLois Curfman McInnes /* wait on receives */ 2897ef1d9bdSSatish Balay while (1) { 2908798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2917ef1d9bdSSatish Balay if (!flg) break; 2928965ea79SLois Curfman McInnes 2937ef1d9bdSSatish Balay for (i=0; i<n;) { 2947ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2957ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2967ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2977ef1d9bdSSatish Balay else ncols = n-i; 2987ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2997ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 3007ef1d9bdSSatish Balay i = j; 3018965ea79SLois Curfman McInnes } 3027ef1d9bdSSatish Balay } 3038798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 3048965ea79SLois Curfman McInnes 30539ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 30639ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 3078965ea79SLois Curfman McInnes 3086d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 30939ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 3108965ea79SLois Curfman McInnes } 3113a40ed3dSBarry Smith PetscFunctionReturn(0); 3128965ea79SLois Curfman McInnes } 3138965ea79SLois Curfman McInnes 3144a2ae208SSatish Balay #undef __FUNCT__ 3154a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 316dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 3178965ea79SLois Curfman McInnes { 318dfbe8321SBarry Smith PetscErrorCode ierr; 31939ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3203a40ed3dSBarry Smith 3213a40ed3dSBarry Smith PetscFunctionBegin; 3223a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 3233a40ed3dSBarry Smith PetscFunctionReturn(0); 3248965ea79SLois Curfman McInnes } 3258965ea79SLois Curfman McInnes 3268965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 3278965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 3288965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 3293501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 3308965ea79SLois Curfman McInnes routine. 3318965ea79SLois Curfman McInnes */ 3324a2ae208SSatish Balay #undef __FUNCT__ 3334a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 334f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 3358965ea79SLois Curfman McInnes { 33639ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3376849ba73SBarry Smith PetscErrorCode ierr; 338f4df32b1SMatthew Knepley PetscInt i,*owners = l->rowners; 33913f74950SBarry Smith PetscInt *nprocs,j,idx,nsends; 34013f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 34113f74950SBarry Smith PetscInt *rvalues,tag = A->tag,count,base,slen,*source; 34213f74950SBarry Smith PetscInt *lens,*lrows,*values; 34313f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 3448965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 3458965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3468965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 34735d8aa7fSBarry Smith PetscTruth found; 3488965ea79SLois Curfman McInnes 3493a40ed3dSBarry Smith PetscFunctionBegin; 3508965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 35113f74950SBarry Smith ierr = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr); 35213f74950SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 35313f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/ 3548965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3558965ea79SLois Curfman McInnes idx = rows[i]; 35635d8aa7fSBarry Smith found = PETSC_FALSE; 3578965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3588965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 359c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3608965ea79SLois Curfman McInnes } 3618965ea79SLois Curfman McInnes } 36229bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3638965ea79SLois Curfman McInnes } 364c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 3658965ea79SLois Curfman McInnes 3668965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 367c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3688965ea79SLois Curfman McInnes 3698965ea79SLois Curfman McInnes /* post receives: */ 37013f74950SBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr); 371b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3728965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 37313f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3748965ea79SLois Curfman McInnes } 3758965ea79SLois Curfman McInnes 3768965ea79SLois Curfman McInnes /* do sends: 3778965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3788965ea79SLois Curfman McInnes the ith processor 3798965ea79SLois Curfman McInnes */ 38013f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr); 381b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 38213f74950SBarry Smith ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr); 3838965ea79SLois Curfman McInnes starts[0] = 0; 384c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3858965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3868965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3878965ea79SLois Curfman McInnes } 3888965ea79SLois Curfman McInnes 3898965ea79SLois Curfman McInnes starts[0] = 0; 390c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3918965ea79SLois Curfman McInnes count = 0; 3928965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 393c1dc657dSBarry Smith if (nprocs[2*i+1]) { 39413f74950SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3958965ea79SLois Curfman McInnes } 3968965ea79SLois Curfman McInnes } 397606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3988965ea79SLois Curfman McInnes 3998965ea79SLois Curfman McInnes base = owners[rank]; 4008965ea79SLois Curfman McInnes 4018965ea79SLois Curfman McInnes /* wait on receives */ 40213f74950SBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr); 4038965ea79SLois Curfman McInnes source = lens + nrecvs; 4048965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 4058965ea79SLois Curfman McInnes while (count) { 406ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4078965ea79SLois Curfman McInnes /* unpack receives into our local space */ 40813f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 4098965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 4108965ea79SLois Curfman McInnes lens[imdex] = n; 4118965ea79SLois Curfman McInnes slen += n; 4128965ea79SLois Curfman McInnes count--; 4138965ea79SLois Curfman McInnes } 414606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 4158965ea79SLois Curfman McInnes 4168965ea79SLois Curfman McInnes /* move the data into the send scatter */ 41713f74950SBarry Smith ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr); 4188965ea79SLois Curfman McInnes count = 0; 4198965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 4208965ea79SLois Curfman McInnes values = rvalues + i*nmax; 4218965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 4228965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 4238965ea79SLois Curfman McInnes } 4248965ea79SLois Curfman McInnes } 425606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 426606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 427606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 428606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 4298965ea79SLois Curfman McInnes 4308965ea79SLois Curfman McInnes /* actually zap the local rows */ 431f4df32b1SMatthew Knepley ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr); 432606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 4338965ea79SLois Curfman McInnes 4348965ea79SLois Curfman McInnes /* wait on sends */ 4358965ea79SLois Curfman McInnes if (nsends) { 436b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 437ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 438606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 4398965ea79SLois Curfman McInnes } 440606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 441606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 4428965ea79SLois Curfman McInnes 4433a40ed3dSBarry Smith PetscFunctionReturn(0); 4448965ea79SLois Curfman McInnes } 4458965ea79SLois Curfman McInnes 4464a2ae208SSatish Balay #undef __FUNCT__ 4474a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 448dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4498965ea79SLois Curfman McInnes { 45039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 451dfbe8321SBarry Smith PetscErrorCode ierr; 452c456f294SBarry Smith 4533a40ed3dSBarry Smith PetscFunctionBegin; 45443a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 45543a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 45644cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4573a40ed3dSBarry Smith PetscFunctionReturn(0); 4588965ea79SLois Curfman McInnes } 4598965ea79SLois Curfman McInnes 4604a2ae208SSatish Balay #undef __FUNCT__ 4614a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 462dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4638965ea79SLois Curfman McInnes { 46439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 465dfbe8321SBarry Smith PetscErrorCode ierr; 466c456f294SBarry Smith 4673a40ed3dSBarry Smith PetscFunctionBegin; 46843a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 46943a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 47044cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4713a40ed3dSBarry Smith PetscFunctionReturn(0); 4728965ea79SLois Curfman McInnes } 4738965ea79SLois Curfman McInnes 4744a2ae208SSatish Balay #undef __FUNCT__ 4754a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 476dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 477096963f5SLois Curfman McInnes { 478096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 479dfbe8321SBarry Smith PetscErrorCode ierr; 48087828ca2SBarry Smith PetscScalar zero = 0.0; 481096963f5SLois Curfman McInnes 4823a40ed3dSBarry Smith PetscFunctionBegin; 4832dcb1b2aSMatthew Knepley ierr = VecSet(yy,zero);CHKERRQ(ierr); 4847c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 485537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 486537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4873a40ed3dSBarry Smith PetscFunctionReturn(0); 488096963f5SLois Curfman McInnes } 489096963f5SLois Curfman McInnes 4904a2ae208SSatish Balay #undef __FUNCT__ 4914a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 492dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 493096963f5SLois Curfman McInnes { 494096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 495dfbe8321SBarry Smith PetscErrorCode ierr; 496096963f5SLois Curfman McInnes 4973a40ed3dSBarry Smith PetscFunctionBegin; 4983501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4997c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 500537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 501537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 5023a40ed3dSBarry Smith PetscFunctionReturn(0); 503096963f5SLois Curfman McInnes } 504096963f5SLois Curfman McInnes 5054a2ae208SSatish Balay #undef __FUNCT__ 5064a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 507dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 5088965ea79SLois Curfman McInnes { 50939ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 510096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 511dfbe8321SBarry Smith PetscErrorCode ierr; 51213f74950SBarry Smith PetscInt len,i,n,m = A->m,radd; 51387828ca2SBarry Smith PetscScalar *x,zero = 0.0; 514ed3cc1f0SBarry Smith 5153a40ed3dSBarry Smith PetscFunctionBegin; 5162dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 5171ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 518096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 519273d9f13SBarry Smith if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 520273d9f13SBarry Smith len = PetscMin(a->A->m,a->A->n); 5217ddc982cSLois Curfman McInnes radd = a->rstart*m; 52244cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 523096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 524096963f5SLois Curfman McInnes } 5251ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 5263a40ed3dSBarry Smith PetscFunctionReturn(0); 5278965ea79SLois Curfman McInnes } 5288965ea79SLois Curfman McInnes 5294a2ae208SSatish Balay #undef __FUNCT__ 5304a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 531dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 5328965ea79SLois Curfman McInnes { 5333501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 534dfbe8321SBarry Smith PetscErrorCode ierr; 535ed3cc1f0SBarry Smith 5363a40ed3dSBarry Smith PetscFunctionBegin; 53794d884c6SBarry Smith 538aa482453SBarry Smith #if defined(PETSC_USE_LOG) 53977431f27SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N); 5408965ea79SLois Curfman McInnes #endif 5418798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 542606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 5433501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 5443501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 5453501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 546622d7880SLois Curfman McInnes if (mdn->factor) { 547606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 548606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 549606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 550606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 551622d7880SLois Curfman McInnes } 552606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 553901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr); 554901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 5553a40ed3dSBarry Smith PetscFunctionReturn(0); 5568965ea79SLois Curfman McInnes } 55739ddd567SLois Curfman McInnes 5584a2ae208SSatish Balay #undef __FUNCT__ 5594a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 5606849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5618965ea79SLois Curfman McInnes { 56239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 563dfbe8321SBarry Smith PetscErrorCode ierr; 5647056b6fcSBarry Smith 5653a40ed3dSBarry Smith PetscFunctionBegin; 56639ddd567SLois Curfman McInnes if (mdn->size == 1) { 56739ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5688965ea79SLois Curfman McInnes } 56929bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5703a40ed3dSBarry Smith PetscFunctionReturn(0); 5718965ea79SLois Curfman McInnes } 5728965ea79SLois Curfman McInnes 5734a2ae208SSatish Balay #undef __FUNCT__ 5744a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 5756849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5768965ea79SLois Curfman McInnes { 57739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 578dfbe8321SBarry Smith PetscErrorCode ierr; 57932dcc486SBarry Smith PetscMPIInt size = mdn->size,rank = mdn->rank; 580b0a32e0cSBarry Smith PetscViewerType vtype; 58132077d6dSBarry Smith PetscTruth iascii,isdraw; 582b0a32e0cSBarry Smith PetscViewer sviewer; 583f3ef73ceSBarry Smith PetscViewerFormat format; 5848965ea79SLois Curfman McInnes 5853a40ed3dSBarry Smith PetscFunctionBegin; 58632077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 587fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 58832077d6dSBarry Smith if (iascii) { 589b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 590b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 591456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 5924e220ebcSLois Curfman McInnes MatInfo info; 593888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 59477431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m, 59577431f27SBarry Smith (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 596b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5973501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5983a40ed3dSBarry Smith PetscFunctionReturn(0); 599fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 6003a40ed3dSBarry Smith PetscFunctionReturn(0); 6018965ea79SLois Curfman McInnes } 602f1af5d2fSBarry Smith } else if (isdraw) { 603b0a32e0cSBarry Smith PetscDraw draw; 604f1af5d2fSBarry Smith PetscTruth isnull; 605f1af5d2fSBarry Smith 606b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 607b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 608f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 609f1af5d2fSBarry Smith } 61077ed5343SBarry Smith 6118965ea79SLois Curfman McInnes if (size == 1) { 61239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 6133a40ed3dSBarry Smith } else { 6148965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 6158965ea79SLois Curfman McInnes Mat A; 61613f74950SBarry Smith PetscInt M = mat->M,N = mat->N,m,row,i,nz; 617ba8c8a56SBarry Smith PetscInt *cols; 618ba8c8a56SBarry Smith PetscScalar *vals; 6198965ea79SLois Curfman McInnes 620f69a0ea3SMatthew Knepley ierr = MatCreate(mat->comm,&A);CHKERRQ(ierr); 6218965ea79SLois Curfman McInnes if (!rank) { 622f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 6233a40ed3dSBarry Smith } else { 624f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 6258965ea79SLois Curfman McInnes } 626be5d1d56SKris Buschelman /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */ 627878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 628878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 62952e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr); 6308965ea79SLois Curfman McInnes 63139ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 63239ddd567SLois Curfman McInnes but it's quick for now */ 63351022da4SBarry Smith A->insertmode = INSERT_VALUES; 634273d9f13SBarry Smith row = mdn->rstart; m = mdn->A->m; 6358965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 636ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 637ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 638ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 63939ddd567SLois Curfman McInnes row++; 6408965ea79SLois Curfman McInnes } 6418965ea79SLois Curfman McInnes 6426d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6436d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 644b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 645b9b97703SBarry Smith if (!rank) { 6466831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 6478965ea79SLois Curfman McInnes } 648b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 649b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6508965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 6518965ea79SLois Curfman McInnes } 6523a40ed3dSBarry Smith PetscFunctionReturn(0); 6538965ea79SLois Curfman McInnes } 6548965ea79SLois Curfman McInnes 6554a2ae208SSatish Balay #undef __FUNCT__ 6564a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 657dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 6588965ea79SLois Curfman McInnes { 659dfbe8321SBarry Smith PetscErrorCode ierr; 66032077d6dSBarry Smith PetscTruth iascii,isbinary,isdraw,issocket; 6618965ea79SLois Curfman McInnes 662433994e6SBarry Smith PetscFunctionBegin; 6630f5bd95cSBarry Smith 66432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 665fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 666b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 667fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 6680f5bd95cSBarry Smith 66932077d6dSBarry Smith if (iascii || issocket || isdraw) { 670f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 6710f5bd95cSBarry Smith } else if (isbinary) { 6723a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6735cd90555SBarry Smith } else { 674958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6758965ea79SLois Curfman McInnes } 6763a40ed3dSBarry Smith PetscFunctionReturn(0); 6778965ea79SLois Curfman McInnes } 6788965ea79SLois Curfman McInnes 6794a2ae208SSatish Balay #undef __FUNCT__ 6804a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 681dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6828965ea79SLois Curfman McInnes { 6833501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 6843501a2bdSLois Curfman McInnes Mat mdn = mat->A; 685dfbe8321SBarry Smith PetscErrorCode ierr; 686329f5518SBarry Smith PetscReal isend[5],irecv[5]; 6878965ea79SLois Curfman McInnes 6883a40ed3dSBarry Smith PetscFunctionBegin; 689273d9f13SBarry Smith info->rows_global = (double)A->M; 690273d9f13SBarry Smith info->columns_global = (double)A->N; 691273d9f13SBarry Smith info->rows_local = (double)A->m; 692273d9f13SBarry Smith info->columns_local = (double)A->N; 6934e220ebcSLois Curfman McInnes info->block_size = 1.0; 6944e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6954e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6964e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6978965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6984e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6994e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 7004e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 7014e220ebcSLois Curfman McInnes info->memory = isend[3]; 7024e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 7038965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 704d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 7054e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7064e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7074e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7084e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7094e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7108965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 711d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 7124e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7134e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7144e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7154e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7164e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7178965ea79SLois Curfman McInnes } 7184e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 7194e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 7204e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 7213a40ed3dSBarry Smith PetscFunctionReturn(0); 7228965ea79SLois Curfman McInnes } 7238965ea79SLois Curfman McInnes 7244a2ae208SSatish Balay #undef __FUNCT__ 7254a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 726dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op) 7278965ea79SLois Curfman McInnes { 72839ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 729dfbe8321SBarry Smith PetscErrorCode ierr; 7308965ea79SLois Curfman McInnes 7313a40ed3dSBarry Smith PetscFunctionBegin; 73212c028f9SKris Buschelman switch (op) { 73312c028f9SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 73412c028f9SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 73512c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 73612c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 73712c028f9SKris Buschelman case MAT_COLUMNS_SORTED: 73812c028f9SKris Buschelman case MAT_COLUMNS_UNSORTED: 739273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 74012c028f9SKris Buschelman break; 74112c028f9SKris Buschelman case MAT_ROW_ORIENTED: 742273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 743273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 74412c028f9SKris Buschelman break; 74512c028f9SKris Buschelman case MAT_ROWS_SORTED: 74612c028f9SKris Buschelman case MAT_ROWS_UNSORTED: 74712c028f9SKris Buschelman case MAT_YES_NEW_DIAGONALS: 74812c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 74963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatSetOption_MPIDense:Option ignored\n"));CHKERRQ(ierr); 75012c028f9SKris Buschelman break; 75112c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 752273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 753273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 75412c028f9SKris Buschelman break; 75512c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 756273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 75712c028f9SKris Buschelman break; 75812c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 75929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 76077e54ba9SKris Buschelman case MAT_SYMMETRIC: 76177e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 7629a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 7639a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 7649a4540c5SBarry Smith case MAT_HERMITIAN: 7659a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 7669a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 7679a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 76877e54ba9SKris Buschelman break; 76912c028f9SKris Buschelman default: 77029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 7713a40ed3dSBarry Smith } 7723a40ed3dSBarry Smith PetscFunctionReturn(0); 7738965ea79SLois Curfman McInnes } 7748965ea79SLois Curfman McInnes 7758965ea79SLois Curfman McInnes 7764a2ae208SSatish Balay #undef __FUNCT__ 7774a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 778dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7795b2fa520SLois Curfman McInnes { 7805b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7815b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 78287828ca2SBarry Smith PetscScalar *l,*r,x,*v; 783dfbe8321SBarry Smith PetscErrorCode ierr; 78413f74950SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 7855b2fa520SLois Curfman McInnes 7865b2fa520SLois Curfman McInnes PetscFunctionBegin; 78772d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7885b2fa520SLois Curfman McInnes if (ll) { 78972d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 79029bbc08cSBarry Smith if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 7911ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7925b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7935b2fa520SLois Curfman McInnes x = l[i]; 7945b2fa520SLois Curfman McInnes v = mat->v + i; 7955b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7965b2fa520SLois Curfman McInnes } 7971ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 798efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 7995b2fa520SLois Curfman McInnes } 8005b2fa520SLois Curfman McInnes if (rr) { 80172d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 80229bbc08cSBarry Smith if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 8035b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 8045b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 8051ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 8065b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 8075b2fa520SLois Curfman McInnes x = r[i]; 8085b2fa520SLois Curfman McInnes v = mat->v + i*m; 8095b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 8105b2fa520SLois Curfman McInnes } 8111ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 812efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8135b2fa520SLois Curfman McInnes } 8145b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 8155b2fa520SLois Curfman McInnes } 8165b2fa520SLois Curfman McInnes 8174a2ae208SSatish Balay #undef __FUNCT__ 8184a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 819dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 820096963f5SLois Curfman McInnes { 8213501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8223501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 823dfbe8321SBarry Smith PetscErrorCode ierr; 82413f74950SBarry Smith PetscInt i,j; 825329f5518SBarry Smith PetscReal sum = 0.0; 82687828ca2SBarry Smith PetscScalar *v = mat->v; 8273501a2bdSLois Curfman McInnes 8283a40ed3dSBarry Smith PetscFunctionBegin; 8293501a2bdSLois Curfman McInnes if (mdn->size == 1) { 830064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 8313501a2bdSLois Curfman McInnes } else { 8323501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 833273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 834aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 835329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 8363501a2bdSLois Curfman McInnes #else 8373501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 8383501a2bdSLois Curfman McInnes #endif 8393501a2bdSLois Curfman McInnes } 840064f8208SBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 841064f8208SBarry Smith *nrm = sqrt(*nrm); 842efee365bSSatish Balay ierr = PetscLogFlops(2*mdn->A->n*mdn->A->m);CHKERRQ(ierr); 8433a40ed3dSBarry Smith } else if (type == NORM_1) { 844329f5518SBarry Smith PetscReal *tmp,*tmp2; 845b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 846273d9f13SBarry Smith tmp2 = tmp + A->N; 847273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 848064f8208SBarry Smith *nrm = 0.0; 8493501a2bdSLois Curfman McInnes v = mat->v; 850273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 851273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 85267e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8533501a2bdSLois Curfman McInnes } 8543501a2bdSLois Curfman McInnes } 855d7d1e502SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 856273d9f13SBarry Smith for (j=0; j<A->N; j++) { 857064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8583501a2bdSLois Curfman McInnes } 859606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 860efee365bSSatish Balay ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr); 8613a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 862329f5518SBarry Smith PetscReal ntemp; 8633501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 864064f8208SBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 8653a40ed3dSBarry Smith } else { 86629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8673501a2bdSLois Curfman McInnes } 8683501a2bdSLois Curfman McInnes } 8693a40ed3dSBarry Smith PetscFunctionReturn(0); 8703501a2bdSLois Curfman McInnes } 8713501a2bdSLois Curfman McInnes 8724a2ae208SSatish Balay #undef __FUNCT__ 8734a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 874dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout) 8753501a2bdSLois Curfman McInnes { 8763501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8773501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8783501a2bdSLois Curfman McInnes Mat B; 87913f74950SBarry Smith PetscInt M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 8806849ba73SBarry Smith PetscErrorCode ierr; 88113f74950SBarry Smith PetscInt j,i; 88287828ca2SBarry Smith PetscScalar *v; 8833501a2bdSLois Curfman McInnes 8843a40ed3dSBarry Smith PetscFunctionBegin; 8857c922b88SBarry Smith if (!matout && M != N) { 88629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 8877056b6fcSBarry Smith } 888f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 889f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr); 890878740d9SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 891878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr); 8923501a2bdSLois Curfman McInnes 893273d9f13SBarry Smith m = a->A->m; n = a->A->n; v = Aloc->v; 89413f74950SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr); 8953501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 8963501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8973501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8983501a2bdSLois Curfman McInnes v += m; 8993501a2bdSLois Curfman McInnes } 900606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 9016d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9026d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9037c922b88SBarry Smith if (matout) { 9043501a2bdSLois Curfman McInnes *matout = B; 9053501a2bdSLois Curfman McInnes } else { 906273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 9073501a2bdSLois Curfman McInnes } 9083a40ed3dSBarry Smith PetscFunctionReturn(0); 909096963f5SLois Curfman McInnes } 910096963f5SLois Curfman McInnes 911d9eff348SSatish Balay #include "petscblaslapack.h" 9124a2ae208SSatish Balay #undef __FUNCT__ 9134a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 914f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha) 91544cd7ae7SLois Curfman McInnes { 91644cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 91744cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 918f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 9194ce68768SBarry Smith PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N; 920efee365bSSatish Balay PetscErrorCode ierr; 92144cd7ae7SLois Curfman McInnes 9223a40ed3dSBarry Smith PetscFunctionBegin; 923f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 924efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 9253a40ed3dSBarry Smith PetscFunctionReturn(0); 92644cd7ae7SLois Curfman McInnes } 92744cd7ae7SLois Curfman McInnes 9286849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 9298965ea79SLois Curfman McInnes 9304a2ae208SSatish Balay #undef __FUNCT__ 9314a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 932dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A) 933273d9f13SBarry Smith { 934dfbe8321SBarry Smith PetscErrorCode ierr; 935273d9f13SBarry Smith 936273d9f13SBarry Smith PetscFunctionBegin; 937273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 938273d9f13SBarry Smith PetscFunctionReturn(0); 939273d9f13SBarry Smith } 940273d9f13SBarry Smith 9418965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 94209dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 94309dc0095SBarry Smith MatGetRow_MPIDense, 94409dc0095SBarry Smith MatRestoreRow_MPIDense, 94509dc0095SBarry Smith MatMult_MPIDense, 94697304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 9477c922b88SBarry Smith MatMultTranspose_MPIDense, 9487c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 9498965ea79SLois Curfman McInnes 0, 95009dc0095SBarry Smith 0, 95109dc0095SBarry Smith 0, 95297304618SKris Buschelman /*10*/ 0, 95309dc0095SBarry Smith 0, 95409dc0095SBarry Smith 0, 95509dc0095SBarry Smith 0, 95609dc0095SBarry Smith MatTranspose_MPIDense, 95797304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 95897304618SKris Buschelman 0, 95909dc0095SBarry Smith MatGetDiagonal_MPIDense, 9605b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 96109dc0095SBarry Smith MatNorm_MPIDense, 96297304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 96309dc0095SBarry Smith MatAssemblyEnd_MPIDense, 96409dc0095SBarry Smith 0, 96509dc0095SBarry Smith MatSetOption_MPIDense, 96609dc0095SBarry Smith MatZeroEntries_MPIDense, 96797304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 96809dc0095SBarry Smith 0, 96909dc0095SBarry Smith 0, 97009dc0095SBarry Smith 0, 97109dc0095SBarry Smith 0, 97297304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 973273d9f13SBarry Smith 0, 97409dc0095SBarry Smith 0, 97509dc0095SBarry Smith MatGetArray_MPIDense, 97609dc0095SBarry Smith MatRestoreArray_MPIDense, 97797304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 97809dc0095SBarry Smith 0, 97909dc0095SBarry Smith 0, 98009dc0095SBarry Smith 0, 98109dc0095SBarry Smith 0, 98297304618SKris Buschelman /*40*/ 0, 9832ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 98409dc0095SBarry Smith 0, 98509dc0095SBarry Smith MatGetValues_MPIDense, 98609dc0095SBarry Smith 0, 98797304618SKris Buschelman /*45*/ 0, 98809dc0095SBarry Smith MatScale_MPIDense, 98909dc0095SBarry Smith 0, 99009dc0095SBarry Smith 0, 99109dc0095SBarry Smith 0, 992521d7252SBarry Smith /*50*/ 0, 99309dc0095SBarry Smith 0, 99409dc0095SBarry Smith 0, 99509dc0095SBarry Smith 0, 99609dc0095SBarry Smith 0, 99797304618SKris Buschelman /*55*/ 0, 99809dc0095SBarry Smith 0, 99909dc0095SBarry Smith 0, 100009dc0095SBarry Smith 0, 100109dc0095SBarry Smith 0, 100297304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 1003b9b97703SBarry Smith MatDestroy_MPIDense, 1004b9b97703SBarry Smith MatView_MPIDense, 100597304618SKris Buschelman MatGetPetscMaps_Petsc, 100697304618SKris Buschelman 0, 100797304618SKris Buschelman /*65*/ 0, 100897304618SKris Buschelman 0, 100997304618SKris Buschelman 0, 101097304618SKris Buschelman 0, 101197304618SKris Buschelman 0, 101297304618SKris Buschelman /*70*/ 0, 101397304618SKris Buschelman 0, 101497304618SKris Buschelman 0, 101597304618SKris Buschelman 0, 101697304618SKris Buschelman 0, 101797304618SKris Buschelman /*75*/ 0, 101897304618SKris Buschelman 0, 101997304618SKris Buschelman 0, 102097304618SKris Buschelman 0, 102197304618SKris Buschelman 0, 102297304618SKris Buschelman /*80*/ 0, 102397304618SKris Buschelman 0, 102497304618SKris Buschelman 0, 102597304618SKris Buschelman 0, 1026865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense, 1027865e5f61SKris Buschelman 0, 1028865e5f61SKris Buschelman 0, 1029865e5f61SKris Buschelman 0, 1030865e5f61SKris Buschelman 0, 1031865e5f61SKris Buschelman 0, 1032865e5f61SKris Buschelman /*90*/ 0, 1033865e5f61SKris Buschelman 0, 1034865e5f61SKris Buschelman 0, 1035865e5f61SKris Buschelman 0, 1036865e5f61SKris Buschelman 0, 1037865e5f61SKris Buschelman /*95*/ 0, 1038865e5f61SKris Buschelman 0, 1039865e5f61SKris Buschelman 0, 1040865e5f61SKris Buschelman 0}; 10418965ea79SLois Curfman McInnes 1042273d9f13SBarry Smith EXTERN_C_BEGIN 10434a2ae208SSatish Balay #undef __FUNCT__ 1044a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 1045be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1046a23d5eceSKris Buschelman { 1047a23d5eceSKris Buschelman Mat_MPIDense *a; 1048dfbe8321SBarry Smith PetscErrorCode ierr; 1049a23d5eceSKris Buschelman 1050a23d5eceSKris Buschelman PetscFunctionBegin; 1051a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1052a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1053a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1054a23d5eceSKris Buschelman 1055a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 1056f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1057f69a0ea3SMatthew Knepley ierr = MatSetSizes(a->A,mat->m,mat->N,mat->m,mat->N);CHKERRQ(ierr); 10585c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 10595c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 106052e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 1061a23d5eceSKris Buschelman PetscFunctionReturn(0); 1062a23d5eceSKris Buschelman } 1063a23d5eceSKris Buschelman EXTERN_C_END 1064a23d5eceSKris Buschelman 10650bad9183SKris Buschelman /*MC 1066fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 10670bad9183SKris Buschelman 10680bad9183SKris Buschelman Options Database Keys: 10690bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 10700bad9183SKris Buschelman 10710bad9183SKris Buschelman Level: beginner 10720bad9183SKris Buschelman 10730bad9183SKris Buschelman .seealso: MatCreateMPIDense 10740bad9183SKris Buschelman M*/ 10750bad9183SKris Buschelman 1076a23d5eceSKris Buschelman EXTERN_C_BEGIN 1077a23d5eceSKris Buschelman #undef __FUNCT__ 10784a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1079be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat) 1080273d9f13SBarry Smith { 1081273d9f13SBarry Smith Mat_MPIDense *a; 1082dfbe8321SBarry Smith PetscErrorCode ierr; 108313f74950SBarry Smith PetscInt i; 1084273d9f13SBarry Smith 1085273d9f13SBarry Smith PetscFunctionBegin; 1086b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1087b0a32e0cSBarry Smith mat->data = (void*)a; 1088273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1089273d9f13SBarry Smith mat->factor = 0; 1090273d9f13SBarry Smith mat->mapping = 0; 1091273d9f13SBarry Smith 1092273d9f13SBarry Smith a->factor = 0; 1093273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1094273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1095273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1096273d9f13SBarry Smith 1097273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1098273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1099273d9f13SBarry Smith a->nvec = mat->n; 1100273d9f13SBarry Smith 1101273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 1102273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 11038a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 11048a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1105273d9f13SBarry Smith 1106273d9f13SBarry Smith /* build local table of row and column ownerships */ 110713f74950SBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 1108273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 110952e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr); 111013f74950SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1111273d9f13SBarry Smith a->rowners[0] = 0; 1112273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1113273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1114273d9f13SBarry Smith } 1115273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1116273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 111713f74950SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1118273d9f13SBarry Smith a->cowners[0] = 0; 1119273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1120273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1121273d9f13SBarry Smith } 1122273d9f13SBarry Smith 1123273d9f13SBarry Smith /* build cache for off array entries formed */ 1124273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1125273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1126273d9f13SBarry Smith 1127273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1128273d9f13SBarry Smith a->lvec = 0; 1129273d9f13SBarry Smith a->Mvctx = 0; 1130273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1131273d9f13SBarry Smith 1132273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1133273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1134273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1135a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1136a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1137a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1138273d9f13SBarry Smith PetscFunctionReturn(0); 1139273d9f13SBarry Smith } 1140273d9f13SBarry Smith EXTERN_C_END 1141273d9f13SBarry Smith 1142209238afSKris Buschelman /*MC 1143002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1144209238afSKris Buschelman 1145209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1146209238afSKris Buschelman and MATMPIDENSE otherwise. 1147209238afSKris Buschelman 1148209238afSKris Buschelman Options Database Keys: 1149209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1150209238afSKris Buschelman 1151209238afSKris Buschelman Level: beginner 1152209238afSKris Buschelman 1153209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1154209238afSKris Buschelman M*/ 1155209238afSKris Buschelman 1156209238afSKris Buschelman EXTERN_C_BEGIN 1157209238afSKris Buschelman #undef __FUNCT__ 1158209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1159be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A) 1160dfbe8321SBarry Smith { 11616849ba73SBarry Smith PetscErrorCode ierr; 116213f74950SBarry Smith PetscMPIInt size; 1163209238afSKris Buschelman 1164209238afSKris Buschelman PetscFunctionBegin; 1165209238afSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1166209238afSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1167209238afSKris Buschelman if (size == 1) { 1168209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1169209238afSKris Buschelman } else { 1170209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1171209238afSKris Buschelman } 1172209238afSKris Buschelman PetscFunctionReturn(0); 1173209238afSKris Buschelman } 1174209238afSKris Buschelman EXTERN_C_END 1175209238afSKris Buschelman 11764a2ae208SSatish Balay #undef __FUNCT__ 11774a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1178273d9f13SBarry Smith /*@C 1179273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1180273d9f13SBarry Smith 1181273d9f13SBarry Smith Not collective 1182273d9f13SBarry Smith 1183273d9f13SBarry Smith Input Parameters: 1184273d9f13SBarry Smith . A - the matrix 1185273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1186273d9f13SBarry Smith to control all matrix memory allocation. 1187273d9f13SBarry Smith 1188273d9f13SBarry Smith Notes: 1189273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1190273d9f13SBarry Smith storage by columns. 1191273d9f13SBarry Smith 1192273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1193273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1194273d9f13SBarry Smith set data=PETSC_NULL. 1195273d9f13SBarry Smith 1196273d9f13SBarry Smith Level: intermediate 1197273d9f13SBarry Smith 1198273d9f13SBarry Smith .keywords: matrix,dense, parallel 1199273d9f13SBarry Smith 1200273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1201273d9f13SBarry Smith @*/ 1202be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1203273d9f13SBarry Smith { 1204dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar *); 1205273d9f13SBarry Smith 1206273d9f13SBarry Smith PetscFunctionBegin; 1207565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1208a23d5eceSKris Buschelman if (f) { 1209a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1210a23d5eceSKris Buschelman } 1211273d9f13SBarry Smith PetscFunctionReturn(0); 1212273d9f13SBarry Smith } 1213273d9f13SBarry Smith 12144a2ae208SSatish Balay #undef __FUNCT__ 12154a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 12168965ea79SLois Curfman McInnes /*@C 121739ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 12188965ea79SLois Curfman McInnes 1219db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1220db81eaa0SLois Curfman McInnes 12218965ea79SLois Curfman McInnes Input Parameters: 1222db81eaa0SLois Curfman McInnes + comm - MPI communicator 12238965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1224db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 12258965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1226db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 12277f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1228dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 12298965ea79SLois Curfman McInnes 12308965ea79SLois Curfman McInnes Output Parameter: 1231477f1c0bSLois Curfman McInnes . A - the matrix 12328965ea79SLois Curfman McInnes 1233b259b22eSLois Curfman McInnes Notes: 123439ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 123539ddd567SLois Curfman McInnes storage by columns. 12368965ea79SLois Curfman McInnes 123718f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 123818f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 12397f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 124018f449edSLois Curfman McInnes 12418965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 12428965ea79SLois Curfman McInnes (possibly both). 12438965ea79SLois Curfman McInnes 1244027ccd11SLois Curfman McInnes Level: intermediate 1245027ccd11SLois Curfman McInnes 124639ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 12478965ea79SLois Curfman McInnes 124839ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 12498965ea79SLois Curfman McInnes @*/ 1250be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 12518965ea79SLois Curfman McInnes { 12526849ba73SBarry Smith PetscErrorCode ierr; 125313f74950SBarry Smith PetscMPIInt size; 12548965ea79SLois Curfman McInnes 12553a40ed3dSBarry Smith PetscFunctionBegin; 1256f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1257f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1258273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1259273d9f13SBarry Smith if (size > 1) { 1260273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1261273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1262273d9f13SBarry Smith } else { 1263273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1264273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 12658c469469SLois Curfman McInnes } 12663a40ed3dSBarry Smith PetscFunctionReturn(0); 12678965ea79SLois Curfman McInnes } 12688965ea79SLois Curfman McInnes 12694a2ae208SSatish Balay #undef __FUNCT__ 12704a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 12716849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 12728965ea79SLois Curfman McInnes { 12738965ea79SLois Curfman McInnes Mat mat; 12743501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1275dfbe8321SBarry Smith PetscErrorCode ierr; 12768965ea79SLois Curfman McInnes 12773a40ed3dSBarry Smith PetscFunctionBegin; 12788965ea79SLois Curfman McInnes *newmat = 0; 1279f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&mat);CHKERRQ(ierr); 1280f69a0ea3SMatthew Knepley ierr = MatSetSizes(mat,A->m,A->n,A->M,A->N);CHKERRQ(ierr); 1281be5d1d56SKris Buschelman ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr); 1282b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1283b0a32e0cSBarry Smith mat->data = (void*)a; 1284549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 12853501a2bdSLois Curfman McInnes mat->factor = A->factor; 1286c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1287273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 12888965ea79SLois Curfman McInnes 12898965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 12908965ea79SLois Curfman McInnes a->rend = oldmat->rend; 12918965ea79SLois Curfman McInnes a->size = oldmat->size; 12928965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1293e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1294b9b97703SBarry Smith a->nvec = oldmat->nvec; 12953782ba37SSatish Balay a->donotstash = oldmat->donotstash; 129613f74950SBarry Smith ierr = PetscMalloc((a->size+1)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 129752e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr); 129813f74950SBarry Smith ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr); 12998798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 13008965ea79SLois Curfman McInnes 1301329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 13025609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 130352e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 13048965ea79SLois Curfman McInnes *newmat = mat; 13053a40ed3dSBarry Smith PetscFunctionReturn(0); 13068965ea79SLois Curfman McInnes } 13078965ea79SLois Curfman McInnes 1308e090d566SSatish Balay #include "petscsys.h" 13098965ea79SLois Curfman McInnes 13104a2ae208SSatish Balay #undef __FUNCT__ 13114a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 1312f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat) 131390ace30eSBarry Smith { 13146849ba73SBarry Smith PetscErrorCode ierr; 131532dcc486SBarry Smith PetscMPIInt rank,size; 131613f74950SBarry Smith PetscInt *rowners,i,m,nz,j; 131787828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 131890ace30eSBarry Smith MPI_Status status; 131990ace30eSBarry Smith 13203a40ed3dSBarry Smith PetscFunctionBegin; 1321d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1322d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 132390ace30eSBarry Smith 132490ace30eSBarry Smith /* determine ownership of all rows */ 132590ace30eSBarry Smith m = M/size + ((M % size) > rank); 132613f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 132713f74950SBarry Smith ierr = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 132890ace30eSBarry Smith rowners[0] = 0; 132990ace30eSBarry Smith for (i=2; i<=size; i++) { 133090ace30eSBarry Smith rowners[i] += rowners[i-1]; 133190ace30eSBarry Smith } 133290ace30eSBarry Smith 1333f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1334f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1335be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1336878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 133790ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 133890ace30eSBarry Smith 133990ace30eSBarry Smith if (!rank) { 134087828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 134190ace30eSBarry Smith 134290ace30eSBarry Smith /* read in my part of the matrix numerical values */ 13430752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 134490ace30eSBarry Smith 134590ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 134690ace30eSBarry Smith vals_ptr = vals; 134790ace30eSBarry Smith for (i=0; i<m; i++) { 134890ace30eSBarry Smith for (j=0; j<N; j++) { 134990ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 135090ace30eSBarry Smith } 135190ace30eSBarry Smith } 135290ace30eSBarry Smith 135390ace30eSBarry Smith /* read in other processors and ship out */ 135490ace30eSBarry Smith for (i=1; i<size; i++) { 135590ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 13560752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1357ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 135890ace30eSBarry Smith } 13593a40ed3dSBarry Smith } else { 136090ace30eSBarry Smith /* receive numeric values */ 136187828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 136290ace30eSBarry Smith 136390ace30eSBarry Smith /* receive message of values*/ 1364ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 136590ace30eSBarry Smith 136690ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 136790ace30eSBarry Smith vals_ptr = vals; 136890ace30eSBarry Smith for (i=0; i<m; i++) { 136990ace30eSBarry Smith for (j=0; j<N; j++) { 137090ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 137190ace30eSBarry Smith } 137290ace30eSBarry Smith } 137390ace30eSBarry Smith } 1374606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1375606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 13766d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13776d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13783a40ed3dSBarry Smith PetscFunctionReturn(0); 137990ace30eSBarry Smith } 138090ace30eSBarry Smith 13814a2ae208SSatish Balay #undef __FUNCT__ 13824a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1383f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat) 13848965ea79SLois Curfman McInnes { 13858965ea79SLois Curfman McInnes Mat A; 138687828ca2SBarry Smith PetscScalar *vals,*svals; 138719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 13888965ea79SLois Curfman McInnes MPI_Status status; 138913f74950SBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz; 139013f74950SBarry Smith PetscInt header[4],*rowlengths = 0,M,N,*cols; 139113f74950SBarry Smith PetscInt *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols; 139213f74950SBarry Smith PetscInt i,nz,j,rstart,rend; 139313f74950SBarry Smith int fd; 13946849ba73SBarry Smith PetscErrorCode ierr; 13958965ea79SLois Curfman McInnes 13963a40ed3dSBarry Smith PetscFunctionBegin; 1397d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1398d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13998965ea79SLois Curfman McInnes if (!rank) { 1400b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 14010752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1402552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 14038965ea79SLois Curfman McInnes } 14048965ea79SLois Curfman McInnes 140513f74950SBarry Smith ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 140690ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 140790ace30eSBarry Smith 140890ace30eSBarry Smith /* 140990ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 141090ace30eSBarry Smith */ 141190ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1412be5d1d56SKris Buschelman ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 14133a40ed3dSBarry Smith PetscFunctionReturn(0); 141490ace30eSBarry Smith } 141590ace30eSBarry Smith 14168965ea79SLois Curfman McInnes /* determine ownership of all rows */ 14178965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 141813f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 1419ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 14208965ea79SLois Curfman McInnes rowners[0] = 0; 14218965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 14228965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 14238965ea79SLois Curfman McInnes } 14248965ea79SLois Curfman McInnes rstart = rowners[rank]; 14258965ea79SLois Curfman McInnes rend = rowners[rank+1]; 14268965ea79SLois Curfman McInnes 14278965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 142813f74950SBarry Smith ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr); 14298965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 14308965ea79SLois Curfman McInnes if (!rank) { 143113f74950SBarry Smith ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 14320752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 143313f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr); 14348965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 14351466f1e1SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 1436606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1437ca161407SBarry Smith } else { 14381466f1e1SBarry Smith ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 14398965ea79SLois Curfman McInnes } 14408965ea79SLois Curfman McInnes 14418965ea79SLois Curfman McInnes if (!rank) { 14428965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 144313f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr); 144413f74950SBarry Smith ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 14458965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14468965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 14478965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 14488965ea79SLois Curfman McInnes } 14498965ea79SLois Curfman McInnes } 1450606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 14518965ea79SLois Curfman McInnes 14528965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 14538965ea79SLois Curfman McInnes maxnz = 0; 14548965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14550452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 14568965ea79SLois Curfman McInnes } 145713f74950SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr); 14588965ea79SLois Curfman McInnes 14598965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 14608965ea79SLois Curfman McInnes nz = procsnz[0]; 146113f74950SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14620752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 14638965ea79SLois Curfman McInnes 14648965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 14658965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14668965ea79SLois Curfman McInnes nz = procsnz[i]; 14670752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 146813f74950SBarry Smith ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 14698965ea79SLois Curfman McInnes } 1470606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 14713a40ed3dSBarry Smith } else { 14728965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 14738965ea79SLois Curfman McInnes nz = 0; 14748965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14758965ea79SLois Curfman McInnes nz += ourlens[i]; 14768965ea79SLois Curfman McInnes } 147713f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14788965ea79SLois Curfman McInnes 14798965ea79SLois Curfman McInnes /* receive message of column indices*/ 148013f74950SBarry Smith ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 148113f74950SBarry Smith ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 148229bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14838965ea79SLois Curfman McInnes } 14848965ea79SLois Curfman McInnes 14858965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 148613f74950SBarry Smith ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr); 14878965ea79SLois Curfman McInnes jj = 0; 14888965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14898965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 14908965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 14918965ea79SLois Curfman McInnes jj++; 14928965ea79SLois Curfman McInnes } 14938965ea79SLois Curfman McInnes } 14948965ea79SLois Curfman McInnes 14958965ea79SLois Curfman McInnes /* create our matrix */ 14968965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14978965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 14988965ea79SLois Curfman McInnes } 1499f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1500f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1501be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1502878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 15038965ea79SLois Curfman McInnes A = *newmat; 15048965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15058965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 15068965ea79SLois Curfman McInnes } 15078965ea79SLois Curfman McInnes 15088965ea79SLois Curfman McInnes if (!rank) { 150987828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 15108965ea79SLois Curfman McInnes 15118965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 15128965ea79SLois Curfman McInnes nz = procsnz[0]; 15130752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 15148965ea79SLois Curfman McInnes 15158965ea79SLois Curfman McInnes /* insert into matrix */ 15168965ea79SLois Curfman McInnes jj = rstart; 15178965ea79SLois Curfman McInnes smycols = mycols; 15188965ea79SLois Curfman McInnes svals = vals; 15198965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15208965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 15218965ea79SLois Curfman McInnes smycols += ourlens[i]; 15228965ea79SLois Curfman McInnes svals += ourlens[i]; 15238965ea79SLois Curfman McInnes jj++; 15248965ea79SLois Curfman McInnes } 15258965ea79SLois Curfman McInnes 15268965ea79SLois Curfman McInnes /* read in other processors and ship out */ 15278965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 15288965ea79SLois Curfman McInnes nz = procsnz[i]; 15290752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1530ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 15318965ea79SLois Curfman McInnes } 1532606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 15333a40ed3dSBarry Smith } else { 15348965ea79SLois Curfman McInnes /* receive numeric values */ 153584d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 15368965ea79SLois Curfman McInnes 15378965ea79SLois Curfman McInnes /* receive message of values*/ 1538ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1539ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 154029bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 15418965ea79SLois Curfman McInnes 15428965ea79SLois Curfman McInnes /* insert into matrix */ 15438965ea79SLois Curfman McInnes jj = rstart; 15448965ea79SLois Curfman McInnes smycols = mycols; 15458965ea79SLois Curfman McInnes svals = vals; 15468965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15478965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 15488965ea79SLois Curfman McInnes smycols += ourlens[i]; 15498965ea79SLois Curfman McInnes svals += ourlens[i]; 15508965ea79SLois Curfman McInnes jj++; 15518965ea79SLois Curfman McInnes } 15528965ea79SLois Curfman McInnes } 1553606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1554606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1555606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1556606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 15578965ea79SLois Curfman McInnes 15586d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15596d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15603a40ed3dSBarry Smith PetscFunctionReturn(0); 15618965ea79SLois Curfman McInnes } 156290ace30eSBarry Smith 156390ace30eSBarry Smith 156490ace30eSBarry Smith 156590ace30eSBarry Smith 156690ace30eSBarry Smith 1567