1be1d678aSKris Buschelman 2ed3cc1f0SBarry Smith /* 3ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 4ed3cc1f0SBarry Smith */ 5ed3cc1f0SBarry Smith 604fea9ffSBarry Smith 7c6db04a5SJed Brown #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 228e6c10adSSatish Balay Level: intermediate 238e6c10adSSatish 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; 29ace3abfcSBarry Smith PetscBool flg; 30ab92ecdeSBarry Smith 31ab92ecdeSBarry Smith PetscFunctionBegin; 32251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr); 332205254eSKarl Rupp if (flg) *B = mat->A; 342205254eSKarl Rupp else *B = A; 35ab92ecdeSBarry Smith PetscFunctionReturn(0); 36ab92ecdeSBarry Smith } 37ab92ecdeSBarry Smith 38ab92ecdeSBarry Smith #undef __FUNCT__ 39ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense" 40ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 41ba8c8a56SBarry Smith { 42ba8c8a56SBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 43ba8c8a56SBarry Smith PetscErrorCode ierr; 44d0f46423SBarry Smith PetscInt lrow,rstart = A->rmap->rstart,rend = A->rmap->rend; 45ba8c8a56SBarry Smith 46ba8c8a56SBarry Smith PetscFunctionBegin; 47e7e72b3dSBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows"); 48ba8c8a56SBarry Smith lrow = row - rstart; 49ba8c8a56SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt**)idx,(const PetscScalar**)v);CHKERRQ(ierr); 50ba8c8a56SBarry Smith PetscFunctionReturn(0); 51ba8c8a56SBarry Smith } 52ba8c8a56SBarry Smith 53ba8c8a56SBarry Smith #undef __FUNCT__ 54ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense" 55ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 56ba8c8a56SBarry Smith { 57ba8c8a56SBarry Smith PetscErrorCode ierr; 58ba8c8a56SBarry Smith 59ba8c8a56SBarry Smith PetscFunctionBegin; 60ba8c8a56SBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 61ba8c8a56SBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 62ba8c8a56SBarry Smith PetscFunctionReturn(0); 63ba8c8a56SBarry Smith } 64ba8c8a56SBarry Smith 650de54da6SSatish Balay EXTERN_C_BEGIN 664a2ae208SSatish Balay #undef __FUNCT__ 674a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 6811bd1e4dSLisandro Dalcin PetscErrorCode MatGetDiagonalBlock_MPIDense(Mat A,Mat *a) 690de54da6SSatish Balay { 700de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 716849ba73SBarry Smith PetscErrorCode ierr; 72d0f46423SBarry Smith PetscInt m = A->rmap->n,rstart = A->rmap->rstart; 7387828ca2SBarry Smith PetscScalar *array; 740de54da6SSatish Balay MPI_Comm comm; 7511bd1e4dSLisandro Dalcin Mat B; 760de54da6SSatish Balay 770de54da6SSatish Balay PetscFunctionBegin; 78e32f2f54SBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only square matrices supported."); 790de54da6SSatish Balay 8011bd1e4dSLisandro Dalcin ierr = PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);CHKERRQ(ierr); 8111bd1e4dSLisandro Dalcin if (!B) { 820de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 8311bd1e4dSLisandro Dalcin ierr = MatCreate(comm,&B);CHKERRQ(ierr); 8411bd1e4dSLisandro Dalcin ierr = MatSetSizes(B,m,m,m,m);CHKERRQ(ierr); 8511bd1e4dSLisandro Dalcin ierr = MatSetType(B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr); 868c778c55SBarry Smith ierr = MatDenseGetArray(mdn->A,&array);CHKERRQ(ierr); 8711bd1e4dSLisandro Dalcin ierr = MatSeqDenseSetPreallocation(B,array+m*rstart);CHKERRQ(ierr); 888c778c55SBarry Smith ierr = MatDenseRestoreArray(mdn->A,&array);CHKERRQ(ierr); 8911bd1e4dSLisandro Dalcin ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9011bd1e4dSLisandro Dalcin ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9111bd1e4dSLisandro Dalcin ierr = PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);CHKERRQ(ierr); 9211bd1e4dSLisandro Dalcin *a = B; 9311bd1e4dSLisandro Dalcin ierr = MatDestroy(&B);CHKERRQ(ierr); 942205254eSKarl Rupp } else *a = B; 950de54da6SSatish Balay PetscFunctionReturn(0); 960de54da6SSatish Balay } 970de54da6SSatish Balay EXTERN_C_END 980de54da6SSatish Balay 994a2ae208SSatish Balay #undef __FUNCT__ 1004a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 10113f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) 1028965ea79SLois Curfman McInnes { 10339b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 104dfbe8321SBarry Smith PetscErrorCode ierr; 105d0f46423SBarry Smith PetscInt i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row; 106ace3abfcSBarry Smith PetscBool roworiented = A->roworiented; 1078965ea79SLois Curfman McInnes 1083a40ed3dSBarry Smith PetscFunctionBegin; 10971fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 1108965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 1115ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 112e32f2f54SBarry Smith if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 1138965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 1148965ea79SLois Curfman McInnes row = idxm[i] - rstart; 11539b7565bSBarry Smith if (roworiented) { 11639b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 1173a40ed3dSBarry Smith } else { 1188965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 1195ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 120e32f2f54SBarry Smith if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 12139b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 12239b7565bSBarry Smith } 1238965ea79SLois Curfman McInnes } 1242205254eSKarl Rupp } else if (!A->donotstash) { 1255080c13bSMatthew G Knepley mat->assembled = PETSC_FALSE; 12639b7565bSBarry Smith if (roworiented) { 127b400d20cSBarry Smith ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n,PETSC_FALSE);CHKERRQ(ierr); 128d36fbae8SSatish Balay } else { 129b400d20cSBarry Smith ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m,PETSC_FALSE);CHKERRQ(ierr); 13039b7565bSBarry Smith } 131b49de8d1SLois Curfman McInnes } 132b49de8d1SLois Curfman McInnes } 1333a40ed3dSBarry Smith PetscFunctionReturn(0); 134b49de8d1SLois Curfman McInnes } 135b49de8d1SLois Curfman McInnes 1364a2ae208SSatish Balay #undef __FUNCT__ 1374a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 13813f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 139b49de8d1SLois Curfman McInnes { 140b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 141dfbe8321SBarry Smith PetscErrorCode ierr; 142d0f46423SBarry Smith PetscInt i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row; 143b49de8d1SLois Curfman McInnes 1443a40ed3dSBarry Smith PetscFunctionBegin; 145b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 146e32f2f54SBarry Smith if (idxm[i] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */ 147e32f2f54SBarry Smith if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 148b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 149b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 150b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 151e32f2f54SBarry Smith if (idxn[j] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */ 152e7e72b3dSBarry Smith if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 153b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 154b49de8d1SLois Curfman McInnes } 155e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported"); 1568965ea79SLois Curfman McInnes } 1573a40ed3dSBarry Smith PetscFunctionReturn(0); 1588965ea79SLois Curfman McInnes } 1598965ea79SLois Curfman McInnes 1604a2ae208SSatish Balay #undef __FUNCT__ 1618c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_MPIDense" 1628c778c55SBarry Smith PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[]) 163ff14e315SSatish Balay { 164ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 165dfbe8321SBarry Smith PetscErrorCode ierr; 166ff14e315SSatish Balay 1673a40ed3dSBarry Smith PetscFunctionBegin; 1688c778c55SBarry Smith ierr = MatDenseGetArray(a->A,array);CHKERRQ(ierr); 1693a40ed3dSBarry Smith PetscFunctionReturn(0); 170ff14e315SSatish Balay } 171ff14e315SSatish Balay 1724a2ae208SSatish Balay #undef __FUNCT__ 1734a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 1744aa3045dSJed Brown static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 175ca3fa75bSLois Curfman McInnes { 176ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 177ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 1786849ba73SBarry Smith PetscErrorCode ierr; 1794aa3045dSJed Brown PetscInt i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols; 1805d0c19d7SBarry Smith const PetscInt *irow,*icol; 18187828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 182ca3fa75bSLois Curfman McInnes Mat newmat; 1834aa3045dSJed Brown IS iscol_local; 184ca3fa75bSLois Curfman McInnes 185ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 1864aa3045dSJed Brown ierr = ISAllGather(iscol,&iscol_local);CHKERRQ(ierr); 187ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 1884aa3045dSJed Brown ierr = ISGetIndices(iscol_local,&icol);CHKERRQ(ierr); 189b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 190b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 1914aa3045dSJed Brown ierr = ISGetSize(iscol,&Ncols);CHKERRQ(ierr); /* global number of columns, size of iscol_local */ 192ca3fa75bSLois Curfman McInnes 193ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1947eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1957eba5e9cSLois Curfman McInnes original matrix! */ 196ca3fa75bSLois Curfman McInnes 197ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1987eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 199ca3fa75bSLois Curfman McInnes 200ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 201ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 202e32f2f54SBarry Smith /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 2037eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 204ca3fa75bSLois Curfman McInnes newmat = *B; 205ca3fa75bSLois Curfman McInnes } else { 206ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 207ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr); 2084aa3045dSJed Brown ierr = MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);CHKERRQ(ierr); 2097adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 2100298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 211ca3fa75bSLois Curfman McInnes } 212ca3fa75bSLois Curfman McInnes 213ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 214ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 215ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense*)newmatd->A->data)->v; 216ca3fa75bSLois Curfman McInnes 2174aa3045dSJed Brown for (i=0; i<Ncols; i++) { 21825a33276SHong Zhang av = v + ((Mat_SeqDense*)mat->A->data)->lda*icol[i]; 219ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 2207eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 221ca3fa75bSLois Curfman McInnes } 222ca3fa75bSLois Curfman McInnes } 223ca3fa75bSLois Curfman McInnes 224ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 225ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 226ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 227ca3fa75bSLois Curfman McInnes 228ca3fa75bSLois Curfman McInnes /* Free work space */ 229ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 2305bdf786aSShri Abhyankar ierr = ISRestoreIndices(iscol_local,&icol);CHKERRQ(ierr); 23132bb1f2dSLisandro Dalcin ierr = ISDestroy(&iscol_local);CHKERRQ(ierr); 232ca3fa75bSLois Curfman McInnes *B = newmat; 233ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 234ca3fa75bSLois Curfman McInnes } 235ca3fa75bSLois Curfman McInnes 2364a2ae208SSatish Balay #undef __FUNCT__ 2378c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_MPIDense" 2388c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 239ff14e315SSatish Balay { 24073a71a0fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 24173a71a0fSBarry Smith PetscErrorCode ierr; 24273a71a0fSBarry Smith 2433a40ed3dSBarry Smith PetscFunctionBegin; 2448c778c55SBarry Smith ierr = MatDenseRestoreArray(a->A,array);CHKERRQ(ierr); 2453a40ed3dSBarry Smith PetscFunctionReturn(0); 246ff14e315SSatish Balay } 247ff14e315SSatish Balay 2484a2ae208SSatish Balay #undef __FUNCT__ 2494a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 250dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2518965ea79SLois Curfman McInnes { 25239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 253ce94432eSBarry Smith MPI_Comm comm; 254dfbe8321SBarry Smith PetscErrorCode ierr; 25513f74950SBarry Smith PetscInt nstash,reallocs; 2568965ea79SLois Curfman McInnes InsertMode addv; 2578965ea79SLois Curfman McInnes 2583a40ed3dSBarry Smith PetscFunctionBegin; 259ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 2608965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 261c3aae356SJed Brown ierr = MPI_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);CHKERRQ(ierr); 262ce94432eSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 263e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2648965ea79SLois Curfman McInnes 265d0f46423SBarry Smith ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr); 2668798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 267ae15b995SBarry Smith ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr); 2683a40ed3dSBarry Smith PetscFunctionReturn(0); 2698965ea79SLois Curfman McInnes } 2708965ea79SLois Curfman McInnes 2714a2ae208SSatish Balay #undef __FUNCT__ 2724a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 273dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2748965ea79SLois Curfman McInnes { 27539ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2766849ba73SBarry Smith PetscErrorCode ierr; 27713f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 27813f74950SBarry Smith PetscMPIInt n; 27987828ca2SBarry Smith PetscScalar *val; 280e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2818965ea79SLois Curfman McInnes 2823a40ed3dSBarry Smith PetscFunctionBegin; 2838965ea79SLois Curfman McInnes /* wait on receives */ 2847ef1d9bdSSatish Balay while (1) { 2858798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2867ef1d9bdSSatish Balay if (!flg) break; 2878965ea79SLois Curfman McInnes 2887ef1d9bdSSatish Balay for (i=0; i<n;) { 2897ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2902205254eSKarl Rupp for (j=i,rstart=row[j]; j<n; j++) { 2912205254eSKarl Rupp if (row[j] != rstart) break; 2922205254eSKarl Rupp } 2937ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2947ef1d9bdSSatish Balay else ncols = n-i; 2957ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2967ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2977ef1d9bdSSatish Balay i = j; 2988965ea79SLois Curfman McInnes } 2997ef1d9bdSSatish Balay } 3008798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 3018965ea79SLois Curfman McInnes 30239ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 30339ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 3048965ea79SLois Curfman McInnes 3056d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 30639ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 3078965ea79SLois Curfman McInnes } 3083a40ed3dSBarry Smith PetscFunctionReturn(0); 3098965ea79SLois Curfman McInnes } 3108965ea79SLois Curfman McInnes 3114a2ae208SSatish Balay #undef __FUNCT__ 3124a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 313dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 3148965ea79SLois Curfman McInnes { 315dfbe8321SBarry Smith PetscErrorCode ierr; 31639ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3173a40ed3dSBarry Smith 3183a40ed3dSBarry Smith PetscFunctionBegin; 3193a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 3203a40ed3dSBarry Smith PetscFunctionReturn(0); 3218965ea79SLois Curfman McInnes } 3228965ea79SLois Curfman McInnes 3238965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 3248965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 3258965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 3263501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 3278965ea79SLois Curfman McInnes routine. 3288965ea79SLois Curfman McInnes */ 3294a2ae208SSatish Balay #undef __FUNCT__ 3304a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 3312b40b63fSBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 3328965ea79SLois Curfman McInnes { 33339ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3346849ba73SBarry Smith PetscErrorCode ierr; 335d0f46423SBarry Smith PetscInt i,*owners = A->rmap->range; 33613f74950SBarry Smith PetscInt *nprocs,j,idx,nsends; 33713f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 3387adad957SLisandro Dalcin PetscInt *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source; 33913f74950SBarry Smith PetscInt *lens,*lrows,*values; 34013f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 341ce94432eSBarry Smith MPI_Comm comm; 3428965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3438965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 344ace3abfcSBarry Smith PetscBool found; 34597b48c8fSBarry Smith const PetscScalar *xx; 34697b48c8fSBarry Smith PetscScalar *bb; 3478965ea79SLois Curfman McInnes 3483a40ed3dSBarry Smith PetscFunctionBegin; 349ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 350ce94432eSBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(comm,PETSC_ERR_SUP,"Only handles square matrices"); 351b9679d65SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership"); 3528965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 35313f74950SBarry Smith ierr = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr); 35413f74950SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 35513f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/ 3568965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3578965ea79SLois Curfman McInnes idx = rows[i]; 35835d8aa7fSBarry Smith found = PETSC_FALSE; 3598965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3608965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 361c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3628965ea79SLois Curfman McInnes } 3638965ea79SLois Curfman McInnes } 364e32f2f54SBarry Smith if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3658965ea79SLois Curfman McInnes } 3662205254eSKarl Rupp nsends = 0; 3672205254eSKarl Rupp for (i=0; i<size; i++) nsends += nprocs[2*i+1]; 3688965ea79SLois Curfman McInnes 3698965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 370c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3718965ea79SLois Curfman McInnes 3728965ea79SLois Curfman McInnes /* post receives: */ 37313f74950SBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr); 374b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3758965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 37613f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3778965ea79SLois Curfman McInnes } 3788965ea79SLois Curfman McInnes 3798965ea79SLois Curfman McInnes /* do sends: 3808965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3818965ea79SLois Curfman McInnes the ith processor 3828965ea79SLois Curfman McInnes */ 38313f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr); 384b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 38513f74950SBarry Smith ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr); 3868965ea79SLois Curfman McInnes 3878965ea79SLois Curfman McInnes starts[0] = 0; 3882205254eSKarl Rupp for (i=1; i<size; i++) starts[i] = starts[i-1] + nprocs[2*i-2]; 3892205254eSKarl Rupp for (i=0; i<N; i++) svalues[starts[owner[i]]++] = rows[i]; 3902205254eSKarl Rupp 3912205254eSKarl Rupp starts[0] = 0; 3922205254eSKarl Rupp for (i=1; i<size+1; i++) starts[i] = starts[i-1] + nprocs[2*i-2]; 3938965ea79SLois Curfman McInnes count = 0; 3948965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 395c1dc657dSBarry Smith if (nprocs[2*i+1]) { 39613f74950SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3978965ea79SLois Curfman McInnes } 3988965ea79SLois Curfman McInnes } 399606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 4008965ea79SLois Curfman McInnes 4018965ea79SLois Curfman McInnes base = owners[rank]; 4028965ea79SLois Curfman McInnes 4038965ea79SLois Curfman McInnes /* wait on receives */ 40474ed9c26SBarry Smith ierr = PetscMalloc2(nrecvs,PetscInt,&lens,nrecvs,PetscInt,&source);CHKERRQ(ierr); 40574ed9c26SBarry Smith count = nrecvs; 40674ed9c26SBarry Smith slen = 0; 4078965ea79SLois Curfman McInnes while (count) { 408ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4098965ea79SLois Curfman McInnes /* unpack receives into our local space */ 41013f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 4112205254eSKarl Rupp 4128965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 4138965ea79SLois Curfman McInnes lens[imdex] = n; 4148965ea79SLois Curfman McInnes slen += n; 4158965ea79SLois Curfman McInnes count--; 4168965ea79SLois Curfman McInnes } 417606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 4188965ea79SLois Curfman McInnes 4198965ea79SLois Curfman McInnes /* move the data into the send scatter */ 42013f74950SBarry Smith ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr); 4218965ea79SLois Curfman McInnes count = 0; 4228965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 4238965ea79SLois Curfman McInnes values = rvalues + i*nmax; 4248965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 4258965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 4268965ea79SLois Curfman McInnes } 4278965ea79SLois Curfman McInnes } 428606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 42974ed9c26SBarry Smith ierr = PetscFree2(lens,source);CHKERRQ(ierr); 430606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 431606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 4328965ea79SLois Curfman McInnes 43397b48c8fSBarry Smith /* fix right hand side if needed */ 43497b48c8fSBarry Smith if (x && b) { 43597b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 43697b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 43797b48c8fSBarry Smith for (i=0; i<slen; i++) { 43897b48c8fSBarry Smith bb[lrows[i]] = diag*xx[lrows[i]]; 43997b48c8fSBarry Smith } 44097b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 44197b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 44297b48c8fSBarry Smith } 44397b48c8fSBarry Smith 4448965ea79SLois Curfman McInnes /* actually zap the local rows */ 445b9679d65SBarry Smith ierr = MatZeroRows(l->A,slen,lrows,0.0,0,0);CHKERRQ(ierr); 446e2eb51b1SBarry Smith if (diag != 0.0) { 447b9679d65SBarry Smith Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data; 448b9679d65SBarry Smith PetscInt m = ll->lda, i; 449b9679d65SBarry Smith 450b9679d65SBarry Smith for (i=0; i<slen; i++) { 451b9679d65SBarry Smith ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag; 452b9679d65SBarry Smith } 453b9679d65SBarry Smith } 454606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 4558965ea79SLois Curfman McInnes 4568965ea79SLois Curfman McInnes /* wait on sends */ 4578965ea79SLois Curfman McInnes if (nsends) { 458b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 459ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 460606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 4618965ea79SLois Curfman McInnes } 462606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 463606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 4643a40ed3dSBarry Smith PetscFunctionReturn(0); 4658965ea79SLois Curfman McInnes } 4668965ea79SLois Curfman McInnes 4674a2ae208SSatish Balay #undef __FUNCT__ 4684a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 469dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4708965ea79SLois Curfman McInnes { 47139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 472dfbe8321SBarry Smith PetscErrorCode ierr; 473c456f294SBarry Smith 4743a40ed3dSBarry Smith PetscFunctionBegin; 475ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 476ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 47744cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4783a40ed3dSBarry Smith PetscFunctionReturn(0); 4798965ea79SLois Curfman McInnes } 4808965ea79SLois Curfman McInnes 4814a2ae208SSatish Balay #undef __FUNCT__ 4824a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 483dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4848965ea79SLois Curfman McInnes { 48539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 486dfbe8321SBarry Smith PetscErrorCode ierr; 487c456f294SBarry Smith 4883a40ed3dSBarry Smith PetscFunctionBegin; 489ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 490ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 49144cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4923a40ed3dSBarry Smith PetscFunctionReturn(0); 4938965ea79SLois Curfman McInnes } 4948965ea79SLois Curfman McInnes 4954a2ae208SSatish Balay #undef __FUNCT__ 4964a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 497dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 498096963f5SLois Curfman McInnes { 499096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 500dfbe8321SBarry Smith PetscErrorCode ierr; 50187828ca2SBarry Smith PetscScalar zero = 0.0; 502096963f5SLois Curfman McInnes 5033a40ed3dSBarry Smith PetscFunctionBegin; 5042dcb1b2aSMatthew Knepley ierr = VecSet(yy,zero);CHKERRQ(ierr); 5057c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 506ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 507ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5083a40ed3dSBarry Smith PetscFunctionReturn(0); 509096963f5SLois Curfman McInnes } 510096963f5SLois Curfman McInnes 5114a2ae208SSatish Balay #undef __FUNCT__ 5124a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 513dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 514096963f5SLois Curfman McInnes { 515096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 516dfbe8321SBarry Smith PetscErrorCode ierr; 517096963f5SLois Curfman McInnes 5183a40ed3dSBarry Smith PetscFunctionBegin; 5193501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 5207c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 521ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 522ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5233a40ed3dSBarry Smith PetscFunctionReturn(0); 524096963f5SLois Curfman McInnes } 525096963f5SLois Curfman McInnes 5264a2ae208SSatish Balay #undef __FUNCT__ 5274a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 528dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 5298965ea79SLois Curfman McInnes { 53039ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 531096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 532dfbe8321SBarry Smith PetscErrorCode ierr; 533d0f46423SBarry Smith PetscInt len,i,n,m = A->rmap->n,radd; 53487828ca2SBarry Smith PetscScalar *x,zero = 0.0; 535ed3cc1f0SBarry Smith 5363a40ed3dSBarry Smith PetscFunctionBegin; 5372dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 5381ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 539096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 540e32f2f54SBarry Smith if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 541d0f46423SBarry Smith len = PetscMin(a->A->rmap->n,a->A->cmap->n); 542d0f46423SBarry Smith radd = A->rmap->rstart*m; 54344cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 544096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 545096963f5SLois Curfman McInnes } 5461ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 5473a40ed3dSBarry Smith PetscFunctionReturn(0); 5488965ea79SLois Curfman McInnes } 5498965ea79SLois Curfman McInnes 5504a2ae208SSatish Balay #undef __FUNCT__ 5514a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 552dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 5538965ea79SLois Curfman McInnes { 5543501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 555dfbe8321SBarry Smith PetscErrorCode ierr; 556ed3cc1f0SBarry Smith 5573a40ed3dSBarry Smith PetscFunctionBegin; 558aa482453SBarry Smith #if defined(PETSC_USE_LOG) 559d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N); 5608965ea79SLois Curfman McInnes #endif 5618798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 5626bf464f9SBarry Smith ierr = MatDestroy(&mdn->A);CHKERRQ(ierr); 5636bf464f9SBarry Smith ierr = VecDestroy(&mdn->lvec);CHKERRQ(ierr); 5646bf464f9SBarry Smith ierr = VecScatterDestroy(&mdn->Mvctx);CHKERRQ(ierr); 56501b82886SBarry Smith 566bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 567dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 568*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C","",NULL);CHKERRQ(ierr); 569*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",NULL);CHKERRQ(ierr); 5700298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","",NULL);CHKERRQ(ierr); 5710298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","",NULL);CHKERRQ(ierr); 5720298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","",NULL);CHKERRQ(ierr); 5733a40ed3dSBarry Smith PetscFunctionReturn(0); 5748965ea79SLois Curfman McInnes } 57539ddd567SLois Curfman McInnes 5764a2ae208SSatish Balay #undef __FUNCT__ 5774a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 5786849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5798965ea79SLois Curfman McInnes { 58039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 581dfbe8321SBarry Smith PetscErrorCode ierr; 582aa05aa95SBarry Smith PetscViewerFormat format; 583aa05aa95SBarry Smith int fd; 584d0f46423SBarry Smith PetscInt header[4],mmax,N = mat->cmap->N,i,j,m,k; 585aa05aa95SBarry Smith PetscMPIInt rank,tag = ((PetscObject)viewer)->tag,size; 586578230a0SSatish Balay PetscScalar *work,*v,*vv; 587aa05aa95SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)mdn->A->data; 5887056b6fcSBarry Smith 5893a40ed3dSBarry Smith PetscFunctionBegin; 59039ddd567SLois Curfman McInnes if (mdn->size == 1) { 59139ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 592aa05aa95SBarry Smith } else { 593aa05aa95SBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 594ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr); 595ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 596aa05aa95SBarry Smith 597aa05aa95SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 598f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 599aa05aa95SBarry Smith 600aa05aa95SBarry Smith if (!rank) { 601aa05aa95SBarry Smith /* store the matrix as a dense matrix */ 6020700a824SBarry Smith header[0] = MAT_FILE_CLASSID; 603d0f46423SBarry Smith header[1] = mat->rmap->N; 604aa05aa95SBarry Smith header[2] = N; 605aa05aa95SBarry Smith header[3] = MATRIX_BINARY_FORMAT_DENSE; 606aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 607aa05aa95SBarry Smith 608aa05aa95SBarry Smith /* get largest work array needed for transposing array */ 609d0f46423SBarry Smith mmax = mat->rmap->n; 610aa05aa95SBarry Smith for (i=1; i<size; i++) { 611d0f46423SBarry Smith mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]); 6128965ea79SLois Curfman McInnes } 613aa05aa95SBarry Smith ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&work);CHKERRQ(ierr); 614aa05aa95SBarry Smith 615aa05aa95SBarry Smith /* write out local array, by rows */ 616d0f46423SBarry Smith m = mat->rmap->n; 617aa05aa95SBarry Smith v = a->v; 618aa05aa95SBarry Smith for (j=0; j<N; j++) { 619aa05aa95SBarry Smith for (i=0; i<m; i++) { 620578230a0SSatish Balay work[j + i*N] = *v++; 621aa05aa95SBarry Smith } 622aa05aa95SBarry Smith } 623aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 624aa05aa95SBarry Smith /* get largest work array to receive messages from other processes, excludes process zero */ 625aa05aa95SBarry Smith mmax = 0; 626aa05aa95SBarry Smith for (i=1; i<size; i++) { 627d0f46423SBarry Smith mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]); 628aa05aa95SBarry Smith } 629578230a0SSatish Balay ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&vv);CHKERRQ(ierr); 630aa05aa95SBarry Smith for (k = 1; k < size; k++) { 631f8009846SMatthew Knepley v = vv; 632d0f46423SBarry Smith m = mat->rmap->range[k+1] - mat->rmap->range[k]; 633ce94432eSBarry Smith ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 634aa05aa95SBarry Smith 635aa05aa95SBarry Smith for (j = 0; j < N; j++) { 636aa05aa95SBarry Smith for (i = 0; i < m; i++) { 637578230a0SSatish Balay work[j + i*N] = *v++; 638aa05aa95SBarry Smith } 639aa05aa95SBarry Smith } 640aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 641aa05aa95SBarry Smith } 642aa05aa95SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 643578230a0SSatish Balay ierr = PetscFree(vv);CHKERRQ(ierr); 644aa05aa95SBarry Smith } else { 645ce94432eSBarry Smith ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 646aa05aa95SBarry Smith } 647eb3f98d2SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE)"); 648aa05aa95SBarry Smith } 6493a40ed3dSBarry Smith PetscFunctionReturn(0); 6508965ea79SLois Curfman McInnes } 6518965ea79SLois Curfman McInnes 6529804daf3SBarry Smith #include <petscdraw.h> 6534a2ae208SSatish Balay #undef __FUNCT__ 6544a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 6556849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 6568965ea79SLois Curfman McInnes { 65739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 658dfbe8321SBarry Smith PetscErrorCode ierr; 65932dcc486SBarry Smith PetscMPIInt size = mdn->size,rank = mdn->rank; 66019fd82e9SBarry Smith PetscViewerType vtype; 661ace3abfcSBarry Smith PetscBool iascii,isdraw; 662b0a32e0cSBarry Smith PetscViewer sviewer; 663f3ef73ceSBarry Smith PetscViewerFormat format; 6648965ea79SLois Curfman McInnes 6653a40ed3dSBarry Smith PetscFunctionBegin; 666251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 667251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 66832077d6dSBarry Smith if (iascii) { 669b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 670b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 671456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 6724e220ebcSLois Curfman McInnes MatInfo info; 673888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 6747b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 6757b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 676b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6777b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 6785d00a290SHong Zhang ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 6793a40ed3dSBarry Smith PetscFunctionReturn(0); 680fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 6813a40ed3dSBarry Smith PetscFunctionReturn(0); 6828965ea79SLois Curfman McInnes } 683f1af5d2fSBarry Smith } else if (isdraw) { 684b0a32e0cSBarry Smith PetscDraw draw; 685ace3abfcSBarry Smith PetscBool isnull; 686f1af5d2fSBarry Smith 687b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 688b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 689f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 690f1af5d2fSBarry Smith } 69177ed5343SBarry Smith 6928965ea79SLois Curfman McInnes if (size == 1) { 69339ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 6943a40ed3dSBarry Smith } else { 6958965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 6968965ea79SLois Curfman McInnes Mat A; 697d0f46423SBarry Smith PetscInt M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz; 698ba8c8a56SBarry Smith PetscInt *cols; 699ba8c8a56SBarry Smith PetscScalar *vals; 7008965ea79SLois Curfman McInnes 701ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr); 7028965ea79SLois Curfman McInnes if (!rank) { 703f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 7043a40ed3dSBarry Smith } else { 705f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 7068965ea79SLois Curfman McInnes } 7077adad957SLisandro Dalcin /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */ 708878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 7090298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(A,NULL);CHKERRQ(ierr); 71052e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr); 7118965ea79SLois Curfman McInnes 71239ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 71339ddd567SLois Curfman McInnes but it's quick for now */ 71451022da4SBarry Smith A->insertmode = INSERT_VALUES; 7152205254eSKarl Rupp 7162205254eSKarl Rupp row = mat->rmap->rstart; 7172205254eSKarl Rupp m = mdn->A->rmap->n; 7188965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 719ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 720ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 721ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 72239ddd567SLois Curfman McInnes row++; 7238965ea79SLois Curfman McInnes } 7248965ea79SLois Curfman McInnes 7256d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7266d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 727b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 728b9b97703SBarry Smith if (!rank) { 7297566de4bSShri Abhyankar ierr = PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr); 7307566de4bSShri Abhyankar /* Set the type name to MATMPIDense so that the correct type can be printed out by PetscObjectPrintClassNamePrefixType() in MatView_SeqDense_ASCII()*/ 7317566de4bSShri Abhyankar PetscStrcpy(((PetscObject)((Mat_MPIDense*)(A->data))->A)->type_name,MATMPIDENSE); 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); 7366bf464f9SBarry Smith 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; 746ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw,issocket; 7478965ea79SLois Curfman McInnes 748433994e6SBarry Smith PetscFunctionBegin; 749251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 750251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 751251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr); 752251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 7530f5bd95cSBarry Smith 75432077d6dSBarry Smith if (iascii || issocket || isdraw) { 755f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 7560f5bd95cSBarry Smith } else if (isbinary) { 7573a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 75811aeaf0aSBarry Smith } 7593a40ed3dSBarry Smith PetscFunctionReturn(0); 7608965ea79SLois Curfman McInnes } 7618965ea79SLois Curfman McInnes 7624a2ae208SSatish Balay #undef __FUNCT__ 7634a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 764dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 7658965ea79SLois Curfman McInnes { 7663501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7673501a2bdSLois Curfman McInnes Mat mdn = mat->A; 768dfbe8321SBarry Smith PetscErrorCode ierr; 769329f5518SBarry Smith PetscReal isend[5],irecv[5]; 7708965ea79SLois Curfman McInnes 7713a40ed3dSBarry Smith PetscFunctionBegin; 7724e220ebcSLois Curfman McInnes info->block_size = 1.0; 7732205254eSKarl Rupp 7744e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 7752205254eSKarl Rupp 7764e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 7774e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 7788965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 7794e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 7804e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 7814e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 7824e220ebcSLois Curfman McInnes info->memory = isend[3]; 7834e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 7848965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 785ce94432eSBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 7862205254eSKarl Rupp 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) { 793ce94432eSBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 7942205254eSKarl Rupp 7954e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7964e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7974e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7984e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7994e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 8008965ea79SLois Curfman McInnes } 8014e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 8024e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 8034e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 8043a40ed3dSBarry Smith PetscFunctionReturn(0); 8058965ea79SLois Curfman McInnes } 8068965ea79SLois Curfman McInnes 8074a2ae208SSatish Balay #undef __FUNCT__ 8084a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 809ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg) 8108965ea79SLois Curfman McInnes { 81139ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 812dfbe8321SBarry Smith PetscErrorCode ierr; 8138965ea79SLois Curfman McInnes 8143a40ed3dSBarry Smith PetscFunctionBegin; 81512c028f9SKris Buschelman switch (op) { 816512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 81712c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 81812c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 8194e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 82012c028f9SKris Buschelman break; 82112c028f9SKris Buschelman case MAT_ROW_ORIENTED: 8224e0d8c25SBarry Smith a->roworiented = flg; 8232205254eSKarl Rupp 8244e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 82512c028f9SKris Buschelman break; 8264e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 82713fa8e87SLisandro Dalcin case MAT_KEEP_NONZERO_PATTERN: 82812c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 829290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 83012c028f9SKris Buschelman break; 83112c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 8324e0d8c25SBarry Smith a->donotstash = flg; 83312c028f9SKris Buschelman break; 83477e54ba9SKris Buschelman case MAT_SYMMETRIC: 83577e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 8369a4540c5SBarry Smith case MAT_HERMITIAN: 8379a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 838600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 839290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 84077e54ba9SKris Buschelman break; 84112c028f9SKris Buschelman default: 842e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 8433a40ed3dSBarry Smith } 8443a40ed3dSBarry Smith PetscFunctionReturn(0); 8458965ea79SLois Curfman McInnes } 8468965ea79SLois Curfman McInnes 8478965ea79SLois Curfman McInnes 8484a2ae208SSatish Balay #undef __FUNCT__ 8494a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 850dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 8515b2fa520SLois Curfman McInnes { 8525b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8535b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 85487828ca2SBarry Smith PetscScalar *l,*r,x,*v; 855dfbe8321SBarry Smith PetscErrorCode ierr; 856d0f46423SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n; 8575b2fa520SLois Curfman McInnes 8585b2fa520SLois Curfman McInnes PetscFunctionBegin; 85972d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 8605b2fa520SLois Curfman McInnes if (ll) { 86172d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 862e32f2f54SBarry Smith if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2); 8631ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 8645b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 8655b2fa520SLois Curfman McInnes x = l[i]; 8665b2fa520SLois Curfman McInnes v = mat->v + i; 8675b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 8685b2fa520SLois Curfman McInnes } 8691ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 870efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8715b2fa520SLois Curfman McInnes } 8725b2fa520SLois Curfman McInnes if (rr) { 873175be7b4SMatthew Knepley ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr); 874e32f2f54SBarry Smith if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3); 875ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 876ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8771ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 8785b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 8795b2fa520SLois Curfman McInnes x = r[i]; 8805b2fa520SLois Curfman McInnes v = mat->v + i*m; 8812205254eSKarl Rupp for (j=0; j<m; j++) (*v++) *= x; 8825b2fa520SLois Curfman McInnes } 8831ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 884efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8855b2fa520SLois Curfman McInnes } 8865b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 8875b2fa520SLois Curfman McInnes } 8885b2fa520SLois Curfman McInnes 8894a2ae208SSatish Balay #undef __FUNCT__ 8904a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 891dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 892096963f5SLois Curfman McInnes { 8933501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8943501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 895dfbe8321SBarry Smith PetscErrorCode ierr; 89613f74950SBarry Smith PetscInt i,j; 897329f5518SBarry Smith PetscReal sum = 0.0; 89887828ca2SBarry Smith PetscScalar *v = mat->v; 8993501a2bdSLois Curfman McInnes 9003a40ed3dSBarry Smith PetscFunctionBegin; 9013501a2bdSLois Curfman McInnes if (mdn->size == 1) { 902064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 9033501a2bdSLois Curfman McInnes } else { 9043501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 905d0f46423SBarry Smith for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) { 906329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 9073501a2bdSLois Curfman McInnes } 908ce94432eSBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 9098f1a2a5eSBarry Smith *nrm = PetscSqrtReal(*nrm); 910dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr); 9113a40ed3dSBarry Smith } else if (type == NORM_1) { 912329f5518SBarry Smith PetscReal *tmp,*tmp2; 91374ed9c26SBarry Smith ierr = PetscMalloc2(A->cmap->N,PetscReal,&tmp,A->cmap->N,PetscReal,&tmp2);CHKERRQ(ierr); 91474ed9c26SBarry Smith ierr = PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr); 91574ed9c26SBarry Smith ierr = PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr); 916064f8208SBarry Smith *nrm = 0.0; 9173501a2bdSLois Curfman McInnes v = mat->v; 918d0f46423SBarry Smith for (j=0; j<mdn->A->cmap->n; j++) { 919d0f46423SBarry Smith for (i=0; i<mdn->A->rmap->n; i++) { 92067e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 9213501a2bdSLois Curfman McInnes } 9223501a2bdSLois Curfman McInnes } 923ce94432eSBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 924d0f46423SBarry Smith for (j=0; j<A->cmap->N; j++) { 925064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 9263501a2bdSLois Curfman McInnes } 92774ed9c26SBarry Smith ierr = PetscFree2(tmp,tmp);CHKERRQ(ierr); 928d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 9293a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 930329f5518SBarry Smith PetscReal ntemp; 9313501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 932ce94432eSBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 933ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm"); 9343501a2bdSLois Curfman McInnes } 9353a40ed3dSBarry Smith PetscFunctionReturn(0); 9363501a2bdSLois Curfman McInnes } 9373501a2bdSLois Curfman McInnes 9384a2ae208SSatish Balay #undef __FUNCT__ 9394a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 940fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout) 9413501a2bdSLois Curfman McInnes { 9423501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 9433501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 9443501a2bdSLois Curfman McInnes Mat B; 945d0f46423SBarry Smith PetscInt M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart; 9466849ba73SBarry Smith PetscErrorCode ierr; 94713f74950SBarry Smith PetscInt j,i; 94887828ca2SBarry Smith PetscScalar *v; 9493501a2bdSLois Curfman McInnes 9503a40ed3dSBarry Smith PetscFunctionBegin; 951ce94432eSBarry Smith if (reuse == MAT_REUSE_MATRIX && A == *matout && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place"); 952fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX || A == *matout) { 953ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 954d0f46423SBarry Smith ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr); 9557adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 9560298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 957fc4dec0aSBarry Smith } else { 958fc4dec0aSBarry Smith B = *matout; 959fc4dec0aSBarry Smith } 9603501a2bdSLois Curfman McInnes 961d0f46423SBarry Smith m = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v; 9621acff37aSSatish Balay ierr = PetscMalloc(m*sizeof(PetscInt),&rwork);CHKERRQ(ierr); 9633501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 9641acff37aSSatish Balay for (j=0; j<n; j++) { 9653501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 9663501a2bdSLois Curfman McInnes v += m; 9673501a2bdSLois Curfman McInnes } 968606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 9696d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9706d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 971815cbec1SBarry Smith if (reuse == MAT_INITIAL_MATRIX || *matout != A) { 9723501a2bdSLois Curfman McInnes *matout = B; 9733501a2bdSLois Curfman McInnes } else { 974eb6b5d47SBarry Smith ierr = MatHeaderMerge(A,B);CHKERRQ(ierr); 9753501a2bdSLois Curfman McInnes } 9763a40ed3dSBarry Smith PetscFunctionReturn(0); 977096963f5SLois Curfman McInnes } 978096963f5SLois Curfman McInnes 97944cd7ae7SLois Curfman McInnes 9806849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*); 981d84338a6SBarry Smith extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar); 9828965ea79SLois Curfman McInnes 9834a2ae208SSatish Balay #undef __FUNCT__ 9844994cf47SJed Brown #define __FUNCT__ "MatSetUp_MPIDense" 9854994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A) 986273d9f13SBarry Smith { 987dfbe8321SBarry Smith PetscErrorCode ierr; 988273d9f13SBarry Smith 989273d9f13SBarry Smith PetscFunctionBegin; 990273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 991273d9f13SBarry Smith PetscFunctionReturn(0); 992273d9f13SBarry Smith } 993273d9f13SBarry Smith 99430716080SHong Zhang #undef __FUNCT__ 995488007eeSBarry Smith #define __FUNCT__ "MatAXPY_MPIDense" 996488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 997488007eeSBarry Smith { 998488007eeSBarry Smith PetscErrorCode ierr; 999488007eeSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data; 1000488007eeSBarry Smith 1001488007eeSBarry Smith PetscFunctionBegin; 1002488007eeSBarry Smith ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr); 1003488007eeSBarry Smith PetscFunctionReturn(0); 1004488007eeSBarry Smith } 1005488007eeSBarry Smith 1006ba337c44SJed Brown #undef __FUNCT__ 1007ba337c44SJed Brown #define __FUNCT__ "MatConjugate_MPIDense" 10087087cfbeSBarry Smith PetscErrorCode MatConjugate_MPIDense(Mat mat) 1009ba337c44SJed Brown { 1010ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)mat->data; 1011ba337c44SJed Brown PetscErrorCode ierr; 1012ba337c44SJed Brown 1013ba337c44SJed Brown PetscFunctionBegin; 1014ba337c44SJed Brown ierr = MatConjugate(a->A);CHKERRQ(ierr); 1015ba337c44SJed Brown PetscFunctionReturn(0); 1016ba337c44SJed Brown } 1017ba337c44SJed Brown 1018ba337c44SJed Brown #undef __FUNCT__ 1019ba337c44SJed Brown #define __FUNCT__ "MatRealPart_MPIDense" 1020ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A) 1021ba337c44SJed Brown { 1022ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1023ba337c44SJed Brown PetscErrorCode ierr; 1024ba337c44SJed Brown 1025ba337c44SJed Brown PetscFunctionBegin; 1026ba337c44SJed Brown ierr = MatRealPart(a->A);CHKERRQ(ierr); 1027ba337c44SJed Brown PetscFunctionReturn(0); 1028ba337c44SJed Brown } 1029ba337c44SJed Brown 1030ba337c44SJed Brown #undef __FUNCT__ 1031ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_MPIDense" 1032ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A) 1033ba337c44SJed Brown { 1034ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1035ba337c44SJed Brown PetscErrorCode ierr; 1036ba337c44SJed Brown 1037ba337c44SJed Brown PetscFunctionBegin; 1038ba337c44SJed Brown ierr = MatImaginaryPart(a->A);CHKERRQ(ierr); 1039ba337c44SJed Brown PetscFunctionReturn(0); 1040ba337c44SJed Brown } 1041ba337c44SJed Brown 10420716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*); 10430716a85fSBarry Smith #undef __FUNCT__ 10440716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense" 10450716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms) 10460716a85fSBarry Smith { 10470716a85fSBarry Smith PetscErrorCode ierr; 10480716a85fSBarry Smith PetscInt i,n; 10490716a85fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*) A->data; 10500716a85fSBarry Smith PetscReal *work; 10510716a85fSBarry Smith 10520716a85fSBarry Smith PetscFunctionBegin; 10530298fd71SBarry Smith ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr); 10540716a85fSBarry Smith ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr); 10550716a85fSBarry Smith ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr); 10560716a85fSBarry Smith if (type == NORM_2) { 10570716a85fSBarry Smith for (i=0; i<n; i++) work[i] *= work[i]; 10580716a85fSBarry Smith } 10590716a85fSBarry Smith if (type == NORM_INFINITY) { 10600716a85fSBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr); 10610716a85fSBarry Smith } else { 10620716a85fSBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr); 10630716a85fSBarry Smith } 10640716a85fSBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 10650716a85fSBarry Smith if (type == NORM_2) { 10668f1a2a5eSBarry Smith for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]); 10670716a85fSBarry Smith } 10680716a85fSBarry Smith PetscFunctionReturn(0); 10690716a85fSBarry Smith } 10700716a85fSBarry Smith 107173a71a0fSBarry Smith #undef __FUNCT__ 107273a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_MPIDense" 107373a71a0fSBarry Smith static PetscErrorCode MatSetRandom_MPIDense(Mat x,PetscRandom rctx) 107473a71a0fSBarry Smith { 107573a71a0fSBarry Smith Mat_MPIDense *d = (Mat_MPIDense*)x->data; 107673a71a0fSBarry Smith PetscErrorCode ierr; 107773a71a0fSBarry Smith PetscScalar *a; 107873a71a0fSBarry Smith PetscInt m,n,i; 107973a71a0fSBarry Smith 108073a71a0fSBarry Smith PetscFunctionBegin; 108173a71a0fSBarry Smith ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr); 10828c778c55SBarry Smith ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr); 108373a71a0fSBarry Smith for (i=0; i<m*n; i++) { 108473a71a0fSBarry Smith ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr); 108573a71a0fSBarry Smith } 10868c778c55SBarry Smith ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr); 108773a71a0fSBarry Smith PetscFunctionReturn(0); 108873a71a0fSBarry Smith } 108973a71a0fSBarry Smith 10908965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 109109dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense, 109209dc0095SBarry Smith MatGetRow_MPIDense, 109309dc0095SBarry Smith MatRestoreRow_MPIDense, 109409dc0095SBarry Smith MatMult_MPIDense, 109597304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 10967c922b88SBarry Smith MatMultTranspose_MPIDense, 10977c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 10988965ea79SLois Curfman McInnes 0, 109909dc0095SBarry Smith 0, 110009dc0095SBarry Smith 0, 110197304618SKris Buschelman /* 10*/ 0, 110209dc0095SBarry Smith 0, 110309dc0095SBarry Smith 0, 110409dc0095SBarry Smith 0, 110509dc0095SBarry Smith MatTranspose_MPIDense, 110697304618SKris Buschelman /* 15*/ MatGetInfo_MPIDense, 11076e4ee0c6SHong Zhang MatEqual_MPIDense, 110809dc0095SBarry Smith MatGetDiagonal_MPIDense, 11095b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 111009dc0095SBarry Smith MatNorm_MPIDense, 111197304618SKris Buschelman /* 20*/ MatAssemblyBegin_MPIDense, 111209dc0095SBarry Smith MatAssemblyEnd_MPIDense, 111309dc0095SBarry Smith MatSetOption_MPIDense, 111409dc0095SBarry Smith MatZeroEntries_MPIDense, 1115d519adbfSMatthew Knepley /* 24*/ MatZeroRows_MPIDense, 1116919b68f7SBarry Smith 0, 111701b82886SBarry Smith 0, 111801b82886SBarry Smith 0, 111901b82886SBarry Smith 0, 11204994cf47SJed Brown /* 29*/ MatSetUp_MPIDense, 1121273d9f13SBarry Smith 0, 112209dc0095SBarry Smith 0, 11238c778c55SBarry Smith 0, 11248c778c55SBarry Smith 0, 1125d519adbfSMatthew Knepley /* 34*/ MatDuplicate_MPIDense, 112609dc0095SBarry Smith 0, 112709dc0095SBarry Smith 0, 112809dc0095SBarry Smith 0, 112909dc0095SBarry Smith 0, 1130d519adbfSMatthew Knepley /* 39*/ MatAXPY_MPIDense, 11312ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 113209dc0095SBarry Smith 0, 113309dc0095SBarry Smith MatGetValues_MPIDense, 113409dc0095SBarry Smith 0, 1135d519adbfSMatthew Knepley /* 44*/ 0, 113609dc0095SBarry Smith MatScale_MPIDense, 113709dc0095SBarry Smith 0, 113809dc0095SBarry Smith 0, 113909dc0095SBarry Smith 0, 114073a71a0fSBarry Smith /* 49*/ MatSetRandom_MPIDense, 114109dc0095SBarry Smith 0, 114209dc0095SBarry Smith 0, 114309dc0095SBarry Smith 0, 114409dc0095SBarry Smith 0, 1145d519adbfSMatthew Knepley /* 54*/ 0, 114609dc0095SBarry Smith 0, 114709dc0095SBarry Smith 0, 114809dc0095SBarry Smith 0, 114909dc0095SBarry Smith 0, 1150d519adbfSMatthew Knepley /* 59*/ MatGetSubMatrix_MPIDense, 1151b9b97703SBarry Smith MatDestroy_MPIDense, 1152b9b97703SBarry Smith MatView_MPIDense, 1153357abbc8SBarry Smith 0, 115497304618SKris Buschelman 0, 1155d519adbfSMatthew Knepley /* 64*/ 0, 115697304618SKris Buschelman 0, 115797304618SKris Buschelman 0, 115897304618SKris Buschelman 0, 115997304618SKris Buschelman 0, 1160d519adbfSMatthew Knepley /* 69*/ 0, 116197304618SKris Buschelman 0, 116297304618SKris Buschelman 0, 116397304618SKris Buschelman 0, 116497304618SKris Buschelman 0, 1165d519adbfSMatthew Knepley /* 74*/ 0, 116697304618SKris Buschelman 0, 116797304618SKris Buschelman 0, 116897304618SKris Buschelman 0, 116997304618SKris Buschelman 0, 1170d519adbfSMatthew Knepley /* 79*/ 0, 117197304618SKris Buschelman 0, 117297304618SKris Buschelman 0, 117397304618SKris Buschelman 0, 11745bba2384SShri Abhyankar /* 83*/ MatLoad_MPIDense, 1175865e5f61SKris Buschelman 0, 1176865e5f61SKris Buschelman 0, 1177865e5f61SKris Buschelman 0, 1178865e5f61SKris Buschelman 0, 1179865e5f61SKris Buschelman 0, 1180d519adbfSMatthew Knepley /* 89*/ 1181865e5f61SKris Buschelman 0, 1182865e5f61SKris Buschelman 0, 1183865e5f61SKris Buschelman 0, 11842fbe02b9SBarry Smith 0, 1185ba337c44SJed Brown 0, 1186d519adbfSMatthew Knepley /* 94*/ 0, 1187865e5f61SKris Buschelman 0, 1188865e5f61SKris Buschelman 0, 1189ba337c44SJed Brown 0, 1190ba337c44SJed Brown 0, 1191ba337c44SJed Brown /* 99*/ 0, 1192ba337c44SJed Brown 0, 1193ba337c44SJed Brown 0, 1194ba337c44SJed Brown MatConjugate_MPIDense, 1195ba337c44SJed Brown 0, 1196ba337c44SJed Brown /*104*/ 0, 1197ba337c44SJed Brown MatRealPart_MPIDense, 1198ba337c44SJed Brown MatImaginaryPart_MPIDense, 119986d161a7SShri Abhyankar 0, 120086d161a7SShri Abhyankar 0, 120186d161a7SShri Abhyankar /*109*/ 0, 120286d161a7SShri Abhyankar 0, 120386d161a7SShri Abhyankar 0, 120486d161a7SShri Abhyankar 0, 120586d161a7SShri Abhyankar 0, 120686d161a7SShri Abhyankar /*114*/ 0, 120786d161a7SShri Abhyankar 0, 120886d161a7SShri Abhyankar 0, 120986d161a7SShri Abhyankar 0, 121086d161a7SShri Abhyankar 0, 121186d161a7SShri Abhyankar /*119*/ 0, 121286d161a7SShri Abhyankar 0, 121386d161a7SShri Abhyankar 0, 12140716a85fSBarry Smith 0, 12150716a85fSBarry Smith 0, 12160716a85fSBarry Smith /*124*/ 0, 12170716a85fSBarry Smith MatGetColumnNorms_MPIDense 1218ba337c44SJed Brown }; 12198965ea79SLois Curfman McInnes 1220273d9f13SBarry Smith EXTERN_C_BEGIN 12214a2ae208SSatish Balay #undef __FUNCT__ 1222a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 12237087cfbeSBarry Smith PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1224a23d5eceSKris Buschelman { 1225a23d5eceSKris Buschelman Mat_MPIDense *a; 1226dfbe8321SBarry Smith PetscErrorCode ierr; 1227a23d5eceSKris Buschelman 1228a23d5eceSKris Buschelman PetscFunctionBegin; 1229a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1230a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1231a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1232a23d5eceSKris Buschelman 1233a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 123434ef9618SShri Abhyankar ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 123534ef9618SShri Abhyankar ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 123634ef9618SShri Abhyankar a->nvec = mat->cmap->n; 123734ef9618SShri Abhyankar 1238f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1239d0f46423SBarry Smith ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr); 12405c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 12415c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 124252e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 1243a23d5eceSKris Buschelman PetscFunctionReturn(0); 1244a23d5eceSKris Buschelman } 1245a23d5eceSKris Buschelman EXTERN_C_END 1246a23d5eceSKris Buschelman 1247a23d5eceSKris Buschelman EXTERN_C_BEGIN 1248a23d5eceSKris Buschelman #undef __FUNCT__ 12494a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 12507087cfbeSBarry Smith PetscErrorCode MatCreate_MPIDense(Mat mat) 1251273d9f13SBarry Smith { 1252273d9f13SBarry Smith Mat_MPIDense *a; 1253dfbe8321SBarry Smith PetscErrorCode ierr; 1254273d9f13SBarry Smith 1255273d9f13SBarry Smith PetscFunctionBegin; 125638f2d2fdSLisandro Dalcin ierr = PetscNewLog(mat,Mat_MPIDense,&a);CHKERRQ(ierr); 1257b0a32e0cSBarry Smith mat->data = (void*)a; 1258273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1259273d9f13SBarry Smith 1260273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1261ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);CHKERRQ(ierr); 1262ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);CHKERRQ(ierr); 1263273d9f13SBarry Smith 1264273d9f13SBarry Smith /* build cache for off array entries formed */ 1265273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 12662205254eSKarl Rupp 1267ce94432eSBarry Smith ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);CHKERRQ(ierr); 1268273d9f13SBarry Smith 1269273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1270273d9f13SBarry Smith a->lvec = 0; 1271273d9f13SBarry Smith a->Mvctx = 0; 1272273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1273273d9f13SBarry Smith 1274*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C","MatDenseGetArray_MPIDense",MatDenseGetArray_MPIDense);CHKERRQ(ierr); 1275*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C","MatDenseRestoreArray_MPIDense",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr); 12768c778c55SBarry Smith 1277*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C","MatGetDiagonalBlock_MPIDense",MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1278*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C","MatMPIDenseSetPreallocation_MPIDense",MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1279*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","MatMatMult_MPIAIJ_MPIDense",MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr); 1280*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","MatMatMultSymbolic_MPIAIJ_MPIDense",MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr); 1281*00de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","MatMatMultNumeric_MPIAIJ_MPIDense",MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr); 128238aed534SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr); 1283273d9f13SBarry Smith PetscFunctionReturn(0); 1284273d9f13SBarry Smith } 1285273d9f13SBarry Smith EXTERN_C_END 1286273d9f13SBarry Smith 1287209238afSKris Buschelman /*MC 1288002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1289209238afSKris Buschelman 1290209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1291209238afSKris Buschelman and MATMPIDENSE otherwise. 1292209238afSKris Buschelman 1293209238afSKris Buschelman Options Database Keys: 1294209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1295209238afSKris Buschelman 1296209238afSKris Buschelman Level: beginner 1297209238afSKris Buschelman 129801b82886SBarry Smith 1299209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1300209238afSKris Buschelman M*/ 1301209238afSKris Buschelman 13024a2ae208SSatish Balay #undef __FUNCT__ 13034a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1304273d9f13SBarry Smith /*@C 1305273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1306273d9f13SBarry Smith 1307273d9f13SBarry Smith Not collective 1308273d9f13SBarry Smith 1309273d9f13SBarry Smith Input Parameters: 1310273d9f13SBarry Smith . A - the matrix 13110298fd71SBarry Smith - data - optional location of matrix data. Set data=NULL for PETSc 1312273d9f13SBarry Smith to control all matrix memory allocation. 1313273d9f13SBarry Smith 1314273d9f13SBarry Smith Notes: 1315273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1316273d9f13SBarry Smith storage by columns. 1317273d9f13SBarry Smith 1318273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1319273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 13200298fd71SBarry Smith set data=NULL. 1321273d9f13SBarry Smith 1322273d9f13SBarry Smith Level: intermediate 1323273d9f13SBarry Smith 1324273d9f13SBarry Smith .keywords: matrix,dense, parallel 1325273d9f13SBarry Smith 1326273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1327273d9f13SBarry Smith @*/ 13287087cfbeSBarry Smith PetscErrorCode MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1329273d9f13SBarry Smith { 13304ac538c5SBarry Smith PetscErrorCode ierr; 1331273d9f13SBarry Smith 1332273d9f13SBarry Smith PetscFunctionBegin; 13334ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(mat,data));CHKERRQ(ierr); 1334273d9f13SBarry Smith PetscFunctionReturn(0); 1335273d9f13SBarry Smith } 1336273d9f13SBarry Smith 13374a2ae208SSatish Balay #undef __FUNCT__ 133869b1f4b7SBarry Smith #define __FUNCT__ "MatCreateDense" 13398965ea79SLois Curfman McInnes /*@C 134069b1f4b7SBarry Smith MatCreateDense - Creates a parallel matrix in dense format. 13418965ea79SLois Curfman McInnes 1342db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1343db81eaa0SLois Curfman McInnes 13448965ea79SLois Curfman McInnes Input Parameters: 1345db81eaa0SLois Curfman McInnes + comm - MPI communicator 13468965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1347db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 13488965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1349db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 13500298fd71SBarry Smith - data - optional location of matrix data. Set data=NULL (NULL_SCALAR for Fortran users) for PETSc 1351dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 13528965ea79SLois Curfman McInnes 13538965ea79SLois Curfman McInnes Output Parameter: 1354477f1c0bSLois Curfman McInnes . A - the matrix 13558965ea79SLois Curfman McInnes 1356b259b22eSLois Curfman McInnes Notes: 135739ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 135839ddd567SLois Curfman McInnes storage by columns. 13598965ea79SLois Curfman McInnes 136018f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 136118f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 13620298fd71SBarry Smith set data=NULL (NULL_SCALAR for Fortran users). 136318f449edSLois Curfman McInnes 13648965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 13658965ea79SLois Curfman McInnes (possibly both). 13668965ea79SLois Curfman McInnes 1367027ccd11SLois Curfman McInnes Level: intermediate 1368027ccd11SLois Curfman McInnes 136939ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 13708965ea79SLois Curfman McInnes 137139ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 13728965ea79SLois Curfman McInnes @*/ 137369b1f4b7SBarry Smith PetscErrorCode MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 13748965ea79SLois Curfman McInnes { 13756849ba73SBarry Smith PetscErrorCode ierr; 137613f74950SBarry Smith PetscMPIInt size; 13778965ea79SLois Curfman McInnes 13783a40ed3dSBarry Smith PetscFunctionBegin; 1379f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1380f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1381273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1382273d9f13SBarry Smith if (size > 1) { 1383273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1384273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1385273d9f13SBarry Smith } else { 1386273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1387273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 13888c469469SLois Curfman McInnes } 13893a40ed3dSBarry Smith PetscFunctionReturn(0); 13908965ea79SLois Curfman McInnes } 13918965ea79SLois Curfman McInnes 13924a2ae208SSatish Balay #undef __FUNCT__ 13934a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 13946849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 13958965ea79SLois Curfman McInnes { 13968965ea79SLois Curfman McInnes Mat mat; 13973501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1398dfbe8321SBarry Smith PetscErrorCode ierr; 13998965ea79SLois Curfman McInnes 14003a40ed3dSBarry Smith PetscFunctionBegin; 14018965ea79SLois Curfman McInnes *newmat = 0; 1402ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&mat);CHKERRQ(ierr); 1403d0f46423SBarry Smith ierr = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 14047adad957SLisandro Dalcin ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr); 1405834f8fabSBarry Smith a = (Mat_MPIDense*)mat->data; 1406e04c1aa4SHong Zhang ierr = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 14075aa7edbeSHong Zhang 1408d5f3da31SBarry Smith mat->factortype = A->factortype; 1409c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1410273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 14118965ea79SLois Curfman McInnes 14128965ea79SLois Curfman McInnes a->size = oldmat->size; 14138965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1414e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1415b9b97703SBarry Smith a->nvec = oldmat->nvec; 14163782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1417e04c1aa4SHong Zhang 14181e1e43feSBarry Smith ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr); 14191e1e43feSBarry Smith ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr); 14208965ea79SLois Curfman McInnes 1421329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 14225609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 142352e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 142401b82886SBarry Smith 14258965ea79SLois Curfman McInnes *newmat = mat; 14263a40ed3dSBarry Smith PetscFunctionReturn(0); 14278965ea79SLois Curfman McInnes } 14288965ea79SLois Curfman McInnes 14294a2ae208SSatish Balay #undef __FUNCT__ 14305bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 14315bba2384SShri Abhyankar PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat,PetscInt sizesset) 143286d161a7SShri Abhyankar { 143386d161a7SShri Abhyankar PetscErrorCode ierr; 143486d161a7SShri Abhyankar PetscMPIInt rank,size; 143586d161a7SShri Abhyankar PetscInt *rowners,i,m,nz,j; 143686d161a7SShri Abhyankar PetscScalar *array,*vals,*vals_ptr; 143786d161a7SShri Abhyankar 143886d161a7SShri Abhyankar PetscFunctionBegin; 143986d161a7SShri Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 144086d161a7SShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 144186d161a7SShri Abhyankar 144286d161a7SShri Abhyankar /* determine ownership of all rows */ 144386d161a7SShri Abhyankar if (newmat->rmap->n < 0) m = M/size + ((M % size) > rank); 144486d161a7SShri Abhyankar else m = newmat->rmap->n; 144586d161a7SShri Abhyankar ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 144686d161a7SShri Abhyankar ierr = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 144786d161a7SShri Abhyankar rowners[0] = 0; 144886d161a7SShri Abhyankar for (i=2; i<=size; i++) { 144986d161a7SShri Abhyankar rowners[i] += rowners[i-1]; 145086d161a7SShri Abhyankar } 145186d161a7SShri Abhyankar 145286d161a7SShri Abhyankar if (!sizesset) { 145386d161a7SShri Abhyankar ierr = MatSetSizes(newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 145486d161a7SShri Abhyankar } 14550298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 14568c778c55SBarry Smith ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr); 145786d161a7SShri Abhyankar 145886d161a7SShri Abhyankar if (!rank) { 145986d161a7SShri Abhyankar ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 146086d161a7SShri Abhyankar 146186d161a7SShri Abhyankar /* read in my part of the matrix numerical values */ 146286d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 146386d161a7SShri Abhyankar 146486d161a7SShri Abhyankar /* insert into matrix-by row (this is why cannot directly read into array */ 146586d161a7SShri Abhyankar vals_ptr = vals; 146686d161a7SShri Abhyankar for (i=0; i<m; i++) { 146786d161a7SShri Abhyankar for (j=0; j<N; j++) { 146886d161a7SShri Abhyankar array[i + j*m] = *vals_ptr++; 146986d161a7SShri Abhyankar } 147086d161a7SShri Abhyankar } 147186d161a7SShri Abhyankar 147286d161a7SShri Abhyankar /* read in other processors and ship out */ 147386d161a7SShri Abhyankar for (i=1; i<size; i++) { 147486d161a7SShri Abhyankar nz = (rowners[i+1] - rowners[i])*N; 147586d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1476a25532f0SBarry Smith ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr); 147786d161a7SShri Abhyankar } 147886d161a7SShri Abhyankar } else { 147986d161a7SShri Abhyankar /* receive numeric values */ 148086d161a7SShri Abhyankar ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 148186d161a7SShri Abhyankar 148286d161a7SShri Abhyankar /* receive message of values*/ 1483a25532f0SBarry Smith ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr); 148486d161a7SShri Abhyankar 148586d161a7SShri Abhyankar /* insert into matrix-by row (this is why cannot directly read into array */ 148686d161a7SShri Abhyankar vals_ptr = vals; 148786d161a7SShri Abhyankar for (i=0; i<m; i++) { 148886d161a7SShri Abhyankar for (j=0; j<N; j++) { 148986d161a7SShri Abhyankar array[i + j*m] = *vals_ptr++; 149086d161a7SShri Abhyankar } 149186d161a7SShri Abhyankar } 149286d161a7SShri Abhyankar } 14938c778c55SBarry Smith ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr); 149486d161a7SShri Abhyankar ierr = PetscFree(rowners);CHKERRQ(ierr); 149586d161a7SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 149686d161a7SShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 149786d161a7SShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 149886d161a7SShri Abhyankar PetscFunctionReturn(0); 149986d161a7SShri Abhyankar } 150086d161a7SShri Abhyankar 150186d161a7SShri Abhyankar #undef __FUNCT__ 15025bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_MPIDense" 1503112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer) 150486d161a7SShri Abhyankar { 150586d161a7SShri Abhyankar PetscScalar *vals,*svals; 1506ce94432eSBarry Smith MPI_Comm comm; 150786d161a7SShri Abhyankar MPI_Status status; 150886d161a7SShri Abhyankar PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz; 150986d161a7SShri Abhyankar PetscInt header[4],*rowlengths = 0,M,N,*cols; 151086d161a7SShri Abhyankar PetscInt *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols; 151186d161a7SShri Abhyankar PetscInt i,nz,j,rstart,rend,sizesset=1,grows,gcols; 151286d161a7SShri Abhyankar int fd; 151386d161a7SShri Abhyankar PetscErrorCode ierr; 151486d161a7SShri Abhyankar 151586d161a7SShri Abhyankar PetscFunctionBegin; 1516ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 151786d161a7SShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 151886d161a7SShri Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 151986d161a7SShri Abhyankar if (!rank) { 152086d161a7SShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 152186d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,(char*)header,4,PETSC_INT);CHKERRQ(ierr); 152286d161a7SShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 152386d161a7SShri Abhyankar } 152486d161a7SShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) sizesset = 0; 152586d161a7SShri Abhyankar 152686d161a7SShri Abhyankar ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 152786d161a7SShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 152886d161a7SShri Abhyankar 152986d161a7SShri Abhyankar /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */ 153086d161a7SShri Abhyankar if (sizesset && newmat->rmap->N < 0) newmat->rmap->N = M; 153186d161a7SShri Abhyankar if (sizesset && newmat->cmap->N < 0) newmat->cmap->N = N; 153286d161a7SShri Abhyankar 153386d161a7SShri Abhyankar /* If global sizes are set, check if they are consistent with that given in the file */ 153486d161a7SShri Abhyankar if (sizesset) { 153586d161a7SShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 153686d161a7SShri Abhyankar } 1537abd38a8fSBarry Smith if (sizesset && newmat->rmap->N != grows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%d) and input matrix has (%d)",M,grows); 1538abd38a8fSBarry Smith if (sizesset && newmat->cmap->N != gcols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%d) and input matrix has (%d)",N,gcols); 153986d161a7SShri Abhyankar 154086d161a7SShri Abhyankar /* 154186d161a7SShri Abhyankar Handle case where matrix is stored on disk as a dense matrix 154286d161a7SShri Abhyankar */ 154386d161a7SShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { 15445bba2384SShri Abhyankar ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat,sizesset);CHKERRQ(ierr); 154586d161a7SShri Abhyankar PetscFunctionReturn(0); 154686d161a7SShri Abhyankar } 154786d161a7SShri Abhyankar 154886d161a7SShri Abhyankar /* determine ownership of all rows */ 15492205254eSKarl Rupp if (newmat->rmap->n < 0) { 15502205254eSKarl Rupp ierr = PetscMPIIntCast(M/size + ((M % size) > rank),&m);CHKERRQ(ierr); 15512205254eSKarl Rupp } else { 15522205254eSKarl Rupp ierr = PetscMPIIntCast(newmat->rmap->n,&m);CHKERRQ(ierr); 15532205254eSKarl Rupp } 155486d161a7SShri Abhyankar ierr = PetscMalloc((size+2)*sizeof(PetscMPIInt),&rowners);CHKERRQ(ierr); 155586d161a7SShri Abhyankar ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 155686d161a7SShri Abhyankar rowners[0] = 0; 155786d161a7SShri Abhyankar for (i=2; i<=size; i++) { 155886d161a7SShri Abhyankar rowners[i] += rowners[i-1]; 155986d161a7SShri Abhyankar } 156086d161a7SShri Abhyankar rstart = rowners[rank]; 156186d161a7SShri Abhyankar rend = rowners[rank+1]; 156286d161a7SShri Abhyankar 156386d161a7SShri Abhyankar /* distribute row lengths to all processors */ 156486d161a7SShri Abhyankar ierr = PetscMalloc2(rend-rstart,PetscInt,&ourlens,rend-rstart,PetscInt,&offlens);CHKERRQ(ierr); 156586d161a7SShri Abhyankar if (!rank) { 156686d161a7SShri Abhyankar ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 156786d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 156886d161a7SShri Abhyankar ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr); 156986d161a7SShri Abhyankar for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 157086d161a7SShri Abhyankar ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 157186d161a7SShri Abhyankar ierr = PetscFree(sndcounts);CHKERRQ(ierr); 157286d161a7SShri Abhyankar } else { 157386d161a7SShri Abhyankar ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 157486d161a7SShri Abhyankar } 157586d161a7SShri Abhyankar 157686d161a7SShri Abhyankar if (!rank) { 157786d161a7SShri Abhyankar /* calculate the number of nonzeros on each processor */ 157886d161a7SShri Abhyankar ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr); 157986d161a7SShri Abhyankar ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 158086d161a7SShri Abhyankar for (i=0; i<size; i++) { 158186d161a7SShri Abhyankar for (j=rowners[i]; j< rowners[i+1]; j++) { 158286d161a7SShri Abhyankar procsnz[i] += rowlengths[j]; 158386d161a7SShri Abhyankar } 158486d161a7SShri Abhyankar } 158586d161a7SShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 158686d161a7SShri Abhyankar 158786d161a7SShri Abhyankar /* determine max buffer needed and allocate it */ 158886d161a7SShri Abhyankar maxnz = 0; 158986d161a7SShri Abhyankar for (i=0; i<size; i++) { 159086d161a7SShri Abhyankar maxnz = PetscMax(maxnz,procsnz[i]); 159186d161a7SShri Abhyankar } 159286d161a7SShri Abhyankar ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr); 159386d161a7SShri Abhyankar 159486d161a7SShri Abhyankar /* read in my part of the matrix column indices */ 159586d161a7SShri Abhyankar nz = procsnz[0]; 159686d161a7SShri Abhyankar ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 159786d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 159886d161a7SShri Abhyankar 159986d161a7SShri Abhyankar /* read in every one elses and ship off */ 160086d161a7SShri Abhyankar for (i=1; i<size; i++) { 160186d161a7SShri Abhyankar nz = procsnz[i]; 160286d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 160386d161a7SShri Abhyankar ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 160486d161a7SShri Abhyankar } 160586d161a7SShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 160686d161a7SShri Abhyankar } else { 160786d161a7SShri Abhyankar /* determine buffer space needed for message */ 160886d161a7SShri Abhyankar nz = 0; 160986d161a7SShri Abhyankar for (i=0; i<m; i++) { 161086d161a7SShri Abhyankar nz += ourlens[i]; 161186d161a7SShri Abhyankar } 161286d161a7SShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 161386d161a7SShri Abhyankar 161486d161a7SShri Abhyankar /* receive message of column indices*/ 161586d161a7SShri Abhyankar ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 161686d161a7SShri Abhyankar ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 161786d161a7SShri Abhyankar if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 161886d161a7SShri Abhyankar } 161986d161a7SShri Abhyankar 162086d161a7SShri Abhyankar /* loop over local rows, determining number of off diagonal entries */ 162186d161a7SShri Abhyankar ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr); 162286d161a7SShri Abhyankar jj = 0; 162386d161a7SShri Abhyankar for (i=0; i<m; i++) { 162486d161a7SShri Abhyankar for (j=0; j<ourlens[i]; j++) { 162586d161a7SShri Abhyankar if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 162686d161a7SShri Abhyankar jj++; 162786d161a7SShri Abhyankar } 162886d161a7SShri Abhyankar } 162986d161a7SShri Abhyankar 163086d161a7SShri Abhyankar /* create our matrix */ 16312205254eSKarl Rupp for (i=0; i<m; i++) ourlens[i] -= offlens[i]; 163286d161a7SShri Abhyankar 163386d161a7SShri Abhyankar if (!sizesset) { 163486d161a7SShri Abhyankar ierr = MatSetSizes(newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 163586d161a7SShri Abhyankar } 16360298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 16372205254eSKarl Rupp for (i=0; i<m; i++) ourlens[i] += offlens[i]; 163886d161a7SShri Abhyankar 163986d161a7SShri Abhyankar if (!rank) { 164086d161a7SShri Abhyankar ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 164186d161a7SShri Abhyankar 164286d161a7SShri Abhyankar /* read in my part of the matrix numerical values */ 164386d161a7SShri Abhyankar nz = procsnz[0]; 164486d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 164586d161a7SShri Abhyankar 164686d161a7SShri Abhyankar /* insert into matrix */ 164786d161a7SShri Abhyankar jj = rstart; 164886d161a7SShri Abhyankar smycols = mycols; 164986d161a7SShri Abhyankar svals = vals; 165086d161a7SShri Abhyankar for (i=0; i<m; i++) { 165186d161a7SShri Abhyankar ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 165286d161a7SShri Abhyankar smycols += ourlens[i]; 165386d161a7SShri Abhyankar svals += ourlens[i]; 165486d161a7SShri Abhyankar jj++; 165586d161a7SShri Abhyankar } 165686d161a7SShri Abhyankar 165786d161a7SShri Abhyankar /* read in other processors and ship out */ 165886d161a7SShri Abhyankar for (i=1; i<size; i++) { 165986d161a7SShri Abhyankar nz = procsnz[i]; 166086d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 166186d161a7SShri Abhyankar ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr); 166286d161a7SShri Abhyankar } 166386d161a7SShri Abhyankar ierr = PetscFree(procsnz);CHKERRQ(ierr); 166486d161a7SShri Abhyankar } else { 166586d161a7SShri Abhyankar /* receive numeric values */ 166686d161a7SShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 166786d161a7SShri Abhyankar 166886d161a7SShri Abhyankar /* receive message of values*/ 166986d161a7SShri Abhyankar ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr); 167086d161a7SShri Abhyankar ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 167186d161a7SShri Abhyankar if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 167286d161a7SShri Abhyankar 167386d161a7SShri Abhyankar /* insert into matrix */ 167486d161a7SShri Abhyankar jj = rstart; 167586d161a7SShri Abhyankar smycols = mycols; 167686d161a7SShri Abhyankar svals = vals; 167786d161a7SShri Abhyankar for (i=0; i<m; i++) { 167886d161a7SShri Abhyankar ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 167986d161a7SShri Abhyankar smycols += ourlens[i]; 168086d161a7SShri Abhyankar svals += ourlens[i]; 168186d161a7SShri Abhyankar jj++; 168286d161a7SShri Abhyankar } 168386d161a7SShri Abhyankar } 168486d161a7SShri Abhyankar ierr = PetscFree2(ourlens,offlens);CHKERRQ(ierr); 168586d161a7SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 168686d161a7SShri Abhyankar ierr = PetscFree(mycols);CHKERRQ(ierr); 168786d161a7SShri Abhyankar ierr = PetscFree(rowners);CHKERRQ(ierr); 168886d161a7SShri Abhyankar 168986d161a7SShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169086d161a7SShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169186d161a7SShri Abhyankar PetscFunctionReturn(0); 169286d161a7SShri Abhyankar } 169386d161a7SShri Abhyankar 169486d161a7SShri Abhyankar #undef __FUNCT__ 16956e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense" 1696ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool *flag) 16976e4ee0c6SHong Zhang { 16986e4ee0c6SHong Zhang Mat_MPIDense *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data; 16996e4ee0c6SHong Zhang Mat a,b; 1700ace3abfcSBarry Smith PetscBool flg; 17016e4ee0c6SHong Zhang PetscErrorCode ierr; 170290ace30eSBarry Smith 17036e4ee0c6SHong Zhang PetscFunctionBegin; 17046e4ee0c6SHong Zhang a = matA->A; 17056e4ee0c6SHong Zhang b = matB->A; 17066e4ee0c6SHong Zhang ierr = MatEqual(a,b,&flg);CHKERRQ(ierr); 1707ce94432eSBarry Smith ierr = MPI_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 17086e4ee0c6SHong Zhang PetscFunctionReturn(0); 17096e4ee0c6SHong Zhang } 171090ace30eSBarry Smith 1711