1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3ed3cc1f0SBarry Smith /* 4ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 5ed3cc1f0SBarry Smith */ 6ed3cc1f0SBarry Smith 789ac3759SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" /*I "petscmat.h" I*/ 88965ea79SLois Curfman McInnes 9ba8c8a56SBarry Smith #undef __FUNCT__ 100ad19415SHong Zhang #define __FUNCT__ "MatLUFactorSymbolic_MPIDense" 110ad19415SHong Zhang PetscErrorCode MatLUFactorSymbolic_MPIDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact) 120ad19415SHong Zhang { 130ad19415SHong Zhang PetscErrorCode ierr; 140ad19415SHong Zhang 150ad19415SHong Zhang PetscFunctionBegin; 166017fc9fSHong Zhang ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 170ad19415SHong Zhang PetscFunctionReturn(0); 180ad19415SHong Zhang } 190ad19415SHong Zhang 200ad19415SHong Zhang #undef __FUNCT__ 210ad19415SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_MPIDense" 226017fc9fSHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_MPIDense(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 230ad19415SHong Zhang { 240ad19415SHong Zhang PetscErrorCode ierr; 250ad19415SHong Zhang 260ad19415SHong Zhang PetscFunctionBegin; 276017fc9fSHong Zhang ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 280ad19415SHong Zhang PetscFunctionReturn(0); 290ad19415SHong Zhang } 300ad19415SHong Zhang 310ad19415SHong Zhang #undef __FUNCT__ 32ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix" 33ab92ecdeSBarry Smith /*@ 34ab92ecdeSBarry Smith 35ab92ecdeSBarry Smith MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential 36ab92ecdeSBarry Smith matrix that represents the operator. For sequential matrices it returns itself. 37ab92ecdeSBarry Smith 38ab92ecdeSBarry Smith Input Parameter: 39ab92ecdeSBarry Smith . A - the Seq or MPI dense matrix 40ab92ecdeSBarry Smith 41ab92ecdeSBarry Smith Output Parameter: 42ab92ecdeSBarry Smith . B - the inner matrix 43ab92ecdeSBarry Smith 448e6c10adSSatish Balay Level: intermediate 458e6c10adSSatish Balay 46ab92ecdeSBarry Smith @*/ 47ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B) 48ab92ecdeSBarry Smith { 49ab92ecdeSBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 50ab92ecdeSBarry Smith PetscErrorCode ierr; 51ab92ecdeSBarry Smith PetscTruth flg; 52ab92ecdeSBarry Smith PetscMPIInt size; 53ab92ecdeSBarry Smith 54ab92ecdeSBarry Smith PetscFunctionBegin; 55ab92ecdeSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr); 56ab92ecdeSBarry Smith if (!flg) { /* this check sucks! */ 57ab92ecdeSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATDENSE,&flg);CHKERRQ(ierr); 58ab92ecdeSBarry Smith if (flg) { 59ab92ecdeSBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 60ab92ecdeSBarry Smith if (size == 1) flg = PETSC_FALSE; 61ab92ecdeSBarry Smith } 62ab92ecdeSBarry Smith } 63ab92ecdeSBarry Smith if (flg) { 64ab92ecdeSBarry Smith *B = mat->A; 65ab92ecdeSBarry Smith } else { 66ab92ecdeSBarry Smith *B = A; 67ab92ecdeSBarry Smith } 68ab92ecdeSBarry Smith PetscFunctionReturn(0); 69ab92ecdeSBarry Smith } 70ab92ecdeSBarry Smith 71ab92ecdeSBarry Smith #undef __FUNCT__ 72ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense" 73ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 74ba8c8a56SBarry Smith { 75ba8c8a56SBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 76ba8c8a56SBarry Smith PetscErrorCode ierr; 77ba8c8a56SBarry Smith PetscInt lrow,rstart = mat->rstart,rend = mat->rend; 78ba8c8a56SBarry Smith 79ba8c8a56SBarry Smith PetscFunctionBegin; 80ba8c8a56SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 81ba8c8a56SBarry Smith lrow = row - rstart; 82ba8c8a56SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr); 83ba8c8a56SBarry Smith PetscFunctionReturn(0); 84ba8c8a56SBarry Smith } 85ba8c8a56SBarry Smith 86ba8c8a56SBarry Smith #undef __FUNCT__ 87ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense" 88ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 89ba8c8a56SBarry Smith { 90ba8c8a56SBarry Smith PetscErrorCode ierr; 91ba8c8a56SBarry Smith 92ba8c8a56SBarry Smith PetscFunctionBegin; 93ba8c8a56SBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 94ba8c8a56SBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 95ba8c8a56SBarry Smith PetscFunctionReturn(0); 96ba8c8a56SBarry Smith } 97ba8c8a56SBarry Smith 980de54da6SSatish Balay EXTERN_C_BEGIN 994a2ae208SSatish Balay #undef __FUNCT__ 1004a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 101be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 1020de54da6SSatish Balay { 1030de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 1046849ba73SBarry Smith PetscErrorCode ierr; 10513f74950SBarry Smith PetscInt m = A->m,rstart = mdn->rstart; 10687828ca2SBarry Smith PetscScalar *array; 1070de54da6SSatish Balay MPI_Comm comm; 1080de54da6SSatish Balay 1090de54da6SSatish Balay PetscFunctionBegin; 110273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 1110de54da6SSatish Balay 1120de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 1130de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 1140de54da6SSatish Balay 1150de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 1160de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 117f69a0ea3SMatthew Knepley ierr = MatCreate(comm,B);CHKERRQ(ierr); 118f69a0ea3SMatthew Knepley ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr); 1195c5985e7SKris Buschelman ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr); 1205c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr); 1210de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 1220de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1230de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1240de54da6SSatish Balay 1250de54da6SSatish Balay *iscopy = PETSC_TRUE; 1260de54da6SSatish Balay PetscFunctionReturn(0); 1270de54da6SSatish Balay } 1280de54da6SSatish Balay EXTERN_C_END 1290de54da6SSatish Balay 1304a2ae208SSatish Balay #undef __FUNCT__ 1314a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 13213f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) 1338965ea79SLois Curfman McInnes { 13439b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 135dfbe8321SBarry Smith PetscErrorCode ierr; 13613f74950SBarry Smith PetscInt i,j,rstart = A->rstart,rend = A->rend,row; 137273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 1388965ea79SLois Curfman McInnes 1393a40ed3dSBarry Smith PetscFunctionBegin; 1408965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 1415ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 142273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 1438965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 1448965ea79SLois Curfman McInnes row = idxm[i] - rstart; 14539b7565bSBarry Smith if (roworiented) { 14639b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 1473a40ed3dSBarry Smith } else { 1488965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 1495ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 150273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 15139b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 15239b7565bSBarry Smith } 1538965ea79SLois Curfman McInnes } 1543a40ed3dSBarry Smith } else { 1553782ba37SSatish Balay if (!A->donotstash) { 15639b7565bSBarry Smith if (roworiented) { 1578798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 158d36fbae8SSatish Balay } else { 1598798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 16039b7565bSBarry Smith } 161b49de8d1SLois Curfman McInnes } 162b49de8d1SLois Curfman McInnes } 1633782ba37SSatish Balay } 1643a40ed3dSBarry Smith PetscFunctionReturn(0); 165b49de8d1SLois Curfman McInnes } 166b49de8d1SLois Curfman McInnes 1674a2ae208SSatish Balay #undef __FUNCT__ 1684a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 16913f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 170b49de8d1SLois Curfman McInnes { 171b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 172dfbe8321SBarry Smith PetscErrorCode ierr; 17313f74950SBarry Smith PetscInt i,j,rstart = mdn->rstart,rend = mdn->rend,row; 174b49de8d1SLois Curfman McInnes 1753a40ed3dSBarry Smith PetscFunctionBegin; 176b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 17729bbc08cSBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 178273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 179b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 180b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 181b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 18229bbc08cSBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 183273d9f13SBarry Smith if (idxn[j] >= mat->N) { 18429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 185a8c6a408SBarry Smith } 186b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 187b49de8d1SLois Curfman McInnes } 188a8c6a408SBarry Smith } else { 18929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 1908965ea79SLois Curfman McInnes } 1918965ea79SLois Curfman McInnes } 1923a40ed3dSBarry Smith PetscFunctionReturn(0); 1938965ea79SLois Curfman McInnes } 1948965ea79SLois Curfman McInnes 1954a2ae208SSatish Balay #undef __FUNCT__ 1964a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 197dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 198ff14e315SSatish Balay { 199ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 200dfbe8321SBarry Smith PetscErrorCode ierr; 201ff14e315SSatish Balay 2023a40ed3dSBarry Smith PetscFunctionBegin; 203ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 2043a40ed3dSBarry Smith PetscFunctionReturn(0); 205ff14e315SSatish Balay } 206ff14e315SSatish Balay 2074a2ae208SSatish Balay #undef __FUNCT__ 2084a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 20913f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 210ca3fa75bSLois Curfman McInnes { 211ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 212ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 2136849ba73SBarry Smith PetscErrorCode ierr; 21413f74950SBarry Smith PetscInt i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 21587828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 216ca3fa75bSLois Curfman McInnes Mat newmat; 217ca3fa75bSLois Curfman McInnes 218ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 219ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 220ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 221b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 222b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 223ca3fa75bSLois Curfman McInnes 224ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 2257eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 2267eba5e9cSLois Curfman McInnes original matrix! */ 227ca3fa75bSLois Curfman McInnes 228ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 2297eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 230ca3fa75bSLois Curfman McInnes 231ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 232ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 23329bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 2347eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 235ca3fa75bSLois Curfman McInnes newmat = *B; 236ca3fa75bSLois Curfman McInnes } else { 237ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 238f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr); 239f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr); 240878740d9SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 241878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 242ca3fa75bSLois Curfman McInnes } 243ca3fa75bSLois Curfman McInnes 244ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 245ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 246ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 247ca3fa75bSLois Curfman McInnes 248ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 249ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 250ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 2517eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 252ca3fa75bSLois Curfman McInnes } 253ca3fa75bSLois Curfman McInnes } 254ca3fa75bSLois Curfman McInnes 255ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 256ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 257ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 258ca3fa75bSLois Curfman McInnes 259ca3fa75bSLois Curfman McInnes /* Free work space */ 260ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 261ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 262ca3fa75bSLois Curfman McInnes *B = newmat; 263ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 264ca3fa75bSLois Curfman McInnes } 265ca3fa75bSLois Curfman McInnes 2664a2ae208SSatish Balay #undef __FUNCT__ 2674a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 268dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 269ff14e315SSatish Balay { 2703a40ed3dSBarry Smith PetscFunctionBegin; 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 272ff14e315SSatish Balay } 273ff14e315SSatish Balay 2744a2ae208SSatish Balay #undef __FUNCT__ 2754a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 276dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2778965ea79SLois Curfman McInnes { 27839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 2798965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 280dfbe8321SBarry Smith PetscErrorCode ierr; 28113f74950SBarry Smith PetscInt nstash,reallocs; 2828965ea79SLois Curfman McInnes InsertMode addv; 2838965ea79SLois Curfman McInnes 2843a40ed3dSBarry Smith PetscFunctionBegin; 2858965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 286ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 2877056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 28829bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 2898965ea79SLois Curfman McInnes } 290e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2918965ea79SLois Curfman McInnes 2928798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 2938798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 294*ae15b995SBarry Smith ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr); 2953a40ed3dSBarry Smith PetscFunctionReturn(0); 2968965ea79SLois Curfman McInnes } 2978965ea79SLois Curfman McInnes 2984a2ae208SSatish Balay #undef __FUNCT__ 2994a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 300dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 3018965ea79SLois Curfman McInnes { 30239ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 3036849ba73SBarry Smith PetscErrorCode ierr; 30413f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 30513f74950SBarry Smith PetscMPIInt n; 30687828ca2SBarry Smith PetscScalar *val; 307e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 3088965ea79SLois Curfman McInnes 3093a40ed3dSBarry Smith PetscFunctionBegin; 3108965ea79SLois Curfman McInnes /* wait on receives */ 3117ef1d9bdSSatish Balay while (1) { 3128798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 3137ef1d9bdSSatish Balay if (!flg) break; 3148965ea79SLois Curfman McInnes 3157ef1d9bdSSatish Balay for (i=0; i<n;) { 3167ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 3177ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 3187ef1d9bdSSatish Balay if (j < n) ncols = j-i; 3197ef1d9bdSSatish Balay else ncols = n-i; 3207ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 3217ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 3227ef1d9bdSSatish Balay i = j; 3238965ea79SLois Curfman McInnes } 3247ef1d9bdSSatish Balay } 3258798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 3268965ea79SLois Curfman McInnes 32739ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 32839ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 3298965ea79SLois Curfman McInnes 3306d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 33139ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 3328965ea79SLois Curfman McInnes } 3333a40ed3dSBarry Smith PetscFunctionReturn(0); 3348965ea79SLois Curfman McInnes } 3358965ea79SLois Curfman McInnes 3364a2ae208SSatish Balay #undef __FUNCT__ 3374a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 338dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 3398965ea79SLois Curfman McInnes { 340dfbe8321SBarry Smith PetscErrorCode ierr; 34139ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3423a40ed3dSBarry Smith 3433a40ed3dSBarry Smith PetscFunctionBegin; 3443a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 3453a40ed3dSBarry Smith PetscFunctionReturn(0); 3468965ea79SLois Curfman McInnes } 3478965ea79SLois Curfman McInnes 3488965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 3498965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 3508965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 3513501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 3528965ea79SLois Curfman McInnes routine. 3538965ea79SLois Curfman McInnes */ 3544a2ae208SSatish Balay #undef __FUNCT__ 3554a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 356f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 3578965ea79SLois Curfman McInnes { 35839ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3596849ba73SBarry Smith PetscErrorCode ierr; 360f4df32b1SMatthew Knepley PetscInt i,*owners = l->rowners; 36113f74950SBarry Smith PetscInt *nprocs,j,idx,nsends; 36213f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 36313f74950SBarry Smith PetscInt *rvalues,tag = A->tag,count,base,slen,*source; 36413f74950SBarry Smith PetscInt *lens,*lrows,*values; 36513f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 3668965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 3678965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3688965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 36935d8aa7fSBarry Smith PetscTruth found; 3708965ea79SLois Curfman McInnes 3713a40ed3dSBarry Smith PetscFunctionBegin; 3728965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 37313f74950SBarry Smith ierr = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr); 37413f74950SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 37513f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/ 3768965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3778965ea79SLois Curfman McInnes idx = rows[i]; 37835d8aa7fSBarry Smith found = PETSC_FALSE; 3798965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3808965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 381c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3828965ea79SLois Curfman McInnes } 3838965ea79SLois Curfman McInnes } 38429bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3858965ea79SLois Curfman McInnes } 386c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 3878965ea79SLois Curfman McInnes 3888965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 389c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3908965ea79SLois Curfman McInnes 3918965ea79SLois Curfman McInnes /* post receives: */ 39213f74950SBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr); 393b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3948965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 39513f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3968965ea79SLois Curfman McInnes } 3978965ea79SLois Curfman McInnes 3988965ea79SLois Curfman McInnes /* do sends: 3998965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 4008965ea79SLois Curfman McInnes the ith processor 4018965ea79SLois Curfman McInnes */ 40213f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr); 403b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 40413f74950SBarry Smith ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr); 4058965ea79SLois Curfman McInnes starts[0] = 0; 406c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 4078965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 4088965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 4098965ea79SLois Curfman McInnes } 4108965ea79SLois Curfman McInnes 4118965ea79SLois Curfman McInnes starts[0] = 0; 412c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 4138965ea79SLois Curfman McInnes count = 0; 4148965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 415c1dc657dSBarry Smith if (nprocs[2*i+1]) { 41613f74950SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 4178965ea79SLois Curfman McInnes } 4188965ea79SLois Curfman McInnes } 419606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 4208965ea79SLois Curfman McInnes 4218965ea79SLois Curfman McInnes base = owners[rank]; 4228965ea79SLois Curfman McInnes 4238965ea79SLois Curfman McInnes /* wait on receives */ 42413f74950SBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr); 4258965ea79SLois Curfman McInnes source = lens + nrecvs; 4268965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 4278965ea79SLois Curfman McInnes while (count) { 428ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4298965ea79SLois Curfman McInnes /* unpack receives into our local space */ 43013f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 4318965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 4328965ea79SLois Curfman McInnes lens[imdex] = n; 4338965ea79SLois Curfman McInnes slen += n; 4348965ea79SLois Curfman McInnes count--; 4358965ea79SLois Curfman McInnes } 436606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 4378965ea79SLois Curfman McInnes 4388965ea79SLois Curfman McInnes /* move the data into the send scatter */ 43913f74950SBarry Smith ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr); 4408965ea79SLois Curfman McInnes count = 0; 4418965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 4428965ea79SLois Curfman McInnes values = rvalues + i*nmax; 4438965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 4448965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 4458965ea79SLois Curfman McInnes } 4468965ea79SLois Curfman McInnes } 447606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 448606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 449606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 450606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 4518965ea79SLois Curfman McInnes 4528965ea79SLois Curfman McInnes /* actually zap the local rows */ 453f4df32b1SMatthew Knepley ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr); 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); 4648965ea79SLois Curfman McInnes 4653a40ed3dSBarry Smith PetscFunctionReturn(0); 4668965ea79SLois Curfman McInnes } 4678965ea79SLois Curfman McInnes 4684a2ae208SSatish Balay #undef __FUNCT__ 4694a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 470dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4718965ea79SLois Curfman McInnes { 47239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 473dfbe8321SBarry Smith PetscErrorCode ierr; 474c456f294SBarry Smith 4753a40ed3dSBarry Smith PetscFunctionBegin; 47643a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 47743a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 47844cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4793a40ed3dSBarry Smith PetscFunctionReturn(0); 4808965ea79SLois Curfman McInnes } 4818965ea79SLois Curfman McInnes 4824a2ae208SSatish Balay #undef __FUNCT__ 4834a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 484dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4858965ea79SLois Curfman McInnes { 48639ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 487dfbe8321SBarry Smith PetscErrorCode ierr; 488c456f294SBarry Smith 4893a40ed3dSBarry Smith PetscFunctionBegin; 49043a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 49143a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 49244cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4933a40ed3dSBarry Smith PetscFunctionReturn(0); 4948965ea79SLois Curfman McInnes } 4958965ea79SLois Curfman McInnes 4964a2ae208SSatish Balay #undef __FUNCT__ 4974a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 498dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 499096963f5SLois Curfman McInnes { 500096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 501dfbe8321SBarry Smith PetscErrorCode ierr; 50287828ca2SBarry Smith PetscScalar zero = 0.0; 503096963f5SLois Curfman McInnes 5043a40ed3dSBarry Smith PetscFunctionBegin; 5052dcb1b2aSMatthew Knepley ierr = VecSet(yy,zero);CHKERRQ(ierr); 5067c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 507537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 508537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 5093a40ed3dSBarry Smith PetscFunctionReturn(0); 510096963f5SLois Curfman McInnes } 511096963f5SLois Curfman McInnes 5124a2ae208SSatish Balay #undef __FUNCT__ 5134a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 514dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 515096963f5SLois Curfman McInnes { 516096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 517dfbe8321SBarry Smith PetscErrorCode ierr; 518096963f5SLois Curfman McInnes 5193a40ed3dSBarry Smith PetscFunctionBegin; 5203501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 5217c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 522537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 523537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 5243a40ed3dSBarry Smith PetscFunctionReturn(0); 525096963f5SLois Curfman McInnes } 526096963f5SLois Curfman McInnes 5274a2ae208SSatish Balay #undef __FUNCT__ 5284a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 529dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 5308965ea79SLois Curfman McInnes { 53139ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 532096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 533dfbe8321SBarry Smith PetscErrorCode ierr; 53413f74950SBarry Smith PetscInt len,i,n,m = A->m,radd; 53587828ca2SBarry Smith PetscScalar *x,zero = 0.0; 536ed3cc1f0SBarry Smith 5373a40ed3dSBarry Smith PetscFunctionBegin; 5382dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 5391ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 540096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 541273d9f13SBarry Smith if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 542273d9f13SBarry Smith len = PetscMin(a->A->m,a->A->n); 5437ddc982cSLois Curfman McInnes radd = a->rstart*m; 54444cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 545096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 546096963f5SLois Curfman McInnes } 5471ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 5483a40ed3dSBarry Smith PetscFunctionReturn(0); 5498965ea79SLois Curfman McInnes } 5508965ea79SLois Curfman McInnes 5514a2ae208SSatish Balay #undef __FUNCT__ 5524a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 553dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 5548965ea79SLois Curfman McInnes { 5553501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 556dfbe8321SBarry Smith PetscErrorCode ierr; 557ed3cc1f0SBarry Smith 5583a40ed3dSBarry Smith PetscFunctionBegin; 55994d884c6SBarry Smith 560aa482453SBarry Smith #if defined(PETSC_USE_LOG) 56177431f27SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N); 5628965ea79SLois Curfman McInnes #endif 5638798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 564606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 5653501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 5663501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 5673501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 568622d7880SLois Curfman McInnes if (mdn->factor) { 569606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 570606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 571606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 572606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 573622d7880SLois Curfman McInnes } 574606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 575901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr); 576901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 5773a40ed3dSBarry Smith PetscFunctionReturn(0); 5788965ea79SLois Curfman McInnes } 57939ddd567SLois Curfman McInnes 5804a2ae208SSatish Balay #undef __FUNCT__ 5814a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 5826849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5838965ea79SLois Curfman McInnes { 58439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 585dfbe8321SBarry Smith PetscErrorCode ierr; 5867056b6fcSBarry Smith 5873a40ed3dSBarry Smith PetscFunctionBegin; 58839ddd567SLois Curfman McInnes if (mdn->size == 1) { 58939ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5908965ea79SLois Curfman McInnes } 59129bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5923a40ed3dSBarry Smith PetscFunctionReturn(0); 5938965ea79SLois Curfman McInnes } 5948965ea79SLois Curfman McInnes 5954a2ae208SSatish Balay #undef __FUNCT__ 5964a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 5976849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5988965ea79SLois Curfman McInnes { 59939ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 600dfbe8321SBarry Smith PetscErrorCode ierr; 60132dcc486SBarry Smith PetscMPIInt size = mdn->size,rank = mdn->rank; 602b0a32e0cSBarry Smith PetscViewerType vtype; 60332077d6dSBarry Smith PetscTruth iascii,isdraw; 604b0a32e0cSBarry Smith PetscViewer sviewer; 605f3ef73ceSBarry Smith PetscViewerFormat format; 6068965ea79SLois Curfman McInnes 6073a40ed3dSBarry Smith PetscFunctionBegin; 60832077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 609fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 61032077d6dSBarry Smith if (iascii) { 611b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 612b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 613456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 6144e220ebcSLois Curfman McInnes MatInfo info; 615888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 61677431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m, 61777431f27SBarry Smith (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 618b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6193501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 6203a40ed3dSBarry Smith PetscFunctionReturn(0); 621fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 6223a40ed3dSBarry Smith PetscFunctionReturn(0); 6238965ea79SLois Curfman McInnes } 624f1af5d2fSBarry Smith } else if (isdraw) { 625b0a32e0cSBarry Smith PetscDraw draw; 626f1af5d2fSBarry Smith PetscTruth isnull; 627f1af5d2fSBarry Smith 628b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 629b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 630f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 631f1af5d2fSBarry Smith } 63277ed5343SBarry Smith 6338965ea79SLois Curfman McInnes if (size == 1) { 63439ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 6353a40ed3dSBarry Smith } else { 6368965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 6378965ea79SLois Curfman McInnes Mat A; 63813f74950SBarry Smith PetscInt M = mat->M,N = mat->N,m,row,i,nz; 639ba8c8a56SBarry Smith PetscInt *cols; 640ba8c8a56SBarry Smith PetscScalar *vals; 6418965ea79SLois Curfman McInnes 642f69a0ea3SMatthew Knepley ierr = MatCreate(mat->comm,&A);CHKERRQ(ierr); 6438965ea79SLois Curfman McInnes if (!rank) { 644f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 6453a40ed3dSBarry Smith } else { 646f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 6478965ea79SLois Curfman McInnes } 648be5d1d56SKris Buschelman /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */ 649878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 650878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 65152e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr); 6528965ea79SLois Curfman McInnes 65339ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 65439ddd567SLois Curfman McInnes but it's quick for now */ 65551022da4SBarry Smith A->insertmode = INSERT_VALUES; 656273d9f13SBarry Smith row = mdn->rstart; m = mdn->A->m; 6578965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 658ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 659ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 660ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 66139ddd567SLois Curfman McInnes row++; 6628965ea79SLois Curfman McInnes } 6638965ea79SLois Curfman McInnes 6646d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6656d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 666b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 667b9b97703SBarry Smith if (!rank) { 6686831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 6698965ea79SLois Curfman McInnes } 670b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 671b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6728965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 6738965ea79SLois Curfman McInnes } 6743a40ed3dSBarry Smith PetscFunctionReturn(0); 6758965ea79SLois Curfman McInnes } 6768965ea79SLois Curfman McInnes 6774a2ae208SSatish Balay #undef __FUNCT__ 6784a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 679dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 6808965ea79SLois Curfman McInnes { 681dfbe8321SBarry Smith PetscErrorCode ierr; 68232077d6dSBarry Smith PetscTruth iascii,isbinary,isdraw,issocket; 6838965ea79SLois Curfman McInnes 684433994e6SBarry Smith PetscFunctionBegin; 6850f5bd95cSBarry Smith 68632077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 687fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 688b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 689fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 6900f5bd95cSBarry Smith 69132077d6dSBarry Smith if (iascii || issocket || isdraw) { 692f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 6930f5bd95cSBarry Smith } else if (isbinary) { 6943a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6955cd90555SBarry Smith } else { 696958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6978965ea79SLois Curfman McInnes } 6983a40ed3dSBarry Smith PetscFunctionReturn(0); 6998965ea79SLois Curfman McInnes } 7008965ea79SLois Curfman McInnes 7014a2ae208SSatish Balay #undef __FUNCT__ 7024a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 703dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 7048965ea79SLois Curfman McInnes { 7053501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7063501a2bdSLois Curfman McInnes Mat mdn = mat->A; 707dfbe8321SBarry Smith PetscErrorCode ierr; 708329f5518SBarry Smith PetscReal isend[5],irecv[5]; 7098965ea79SLois Curfman McInnes 7103a40ed3dSBarry Smith PetscFunctionBegin; 711273d9f13SBarry Smith info->rows_global = (double)A->M; 712273d9f13SBarry Smith info->columns_global = (double)A->N; 713273d9f13SBarry Smith info->rows_local = (double)A->m; 714273d9f13SBarry Smith info->columns_local = (double)A->N; 7154e220ebcSLois Curfman McInnes info->block_size = 1.0; 7164e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 7174e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 7184e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 7198965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 7204e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 7214e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 7224e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 7234e220ebcSLois Curfman McInnes info->memory = isend[3]; 7244e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 7258965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 726d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 7274e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7284e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7294e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7304e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7314e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7328965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 733d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 7344e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7354e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7364e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7374e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7384e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7398965ea79SLois Curfman McInnes } 7404e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 7414e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 7424e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 7433a40ed3dSBarry Smith PetscFunctionReturn(0); 7448965ea79SLois Curfman McInnes } 7458965ea79SLois Curfman McInnes 7464a2ae208SSatish Balay #undef __FUNCT__ 7474a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 748dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op) 7498965ea79SLois Curfman McInnes { 75039ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 751dfbe8321SBarry Smith PetscErrorCode ierr; 7528965ea79SLois Curfman McInnes 7533a40ed3dSBarry Smith PetscFunctionBegin; 75412c028f9SKris Buschelman switch (op) { 75512c028f9SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 75612c028f9SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 75712c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 75812c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 75912c028f9SKris Buschelman case MAT_COLUMNS_SORTED: 76012c028f9SKris Buschelman case MAT_COLUMNS_UNSORTED: 761273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 76212c028f9SKris Buschelman break; 76312c028f9SKris Buschelman case MAT_ROW_ORIENTED: 764273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 765273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 76612c028f9SKris Buschelman break; 76712c028f9SKris Buschelman case MAT_ROWS_SORTED: 76812c028f9SKris Buschelman case MAT_ROWS_UNSORTED: 76912c028f9SKris Buschelman case MAT_YES_NEW_DIAGONALS: 77012c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 771*ae15b995SBarry Smith ierr = PetscInfo(A,"Option ignored\n");CHKERRQ(ierr); 77212c028f9SKris Buschelman break; 77312c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 774273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 775273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 77612c028f9SKris Buschelman break; 77712c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 778273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 77912c028f9SKris Buschelman break; 78012c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 78129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 78277e54ba9SKris Buschelman case MAT_SYMMETRIC: 78377e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 7849a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 7859a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 7869a4540c5SBarry Smith case MAT_HERMITIAN: 7879a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 7889a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 7899a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 79077e54ba9SKris Buschelman break; 79112c028f9SKris Buschelman default: 79229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 7933a40ed3dSBarry Smith } 7943a40ed3dSBarry Smith PetscFunctionReturn(0); 7958965ea79SLois Curfman McInnes } 7968965ea79SLois Curfman McInnes 7978965ea79SLois Curfman McInnes 7984a2ae208SSatish Balay #undef __FUNCT__ 7994a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 800dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 8015b2fa520SLois Curfman McInnes { 8025b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8035b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 80487828ca2SBarry Smith PetscScalar *l,*r,x,*v; 805dfbe8321SBarry Smith PetscErrorCode ierr; 80613f74950SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 8075b2fa520SLois Curfman McInnes 8085b2fa520SLois Curfman McInnes PetscFunctionBegin; 80972d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 8105b2fa520SLois Curfman McInnes if (ll) { 81172d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 812175be7b4SMatthew Knepley if (s2a != s2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2); 8131ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 8145b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 8155b2fa520SLois Curfman McInnes x = l[i]; 8165b2fa520SLois Curfman McInnes v = mat->v + i; 8175b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 8185b2fa520SLois Curfman McInnes } 8191ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 820efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8215b2fa520SLois Curfman McInnes } 8225b2fa520SLois Curfman McInnes if (rr) { 823175be7b4SMatthew Knepley ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr); 824175be7b4SMatthew Knepley if (s3a != s3) SETERRQ2(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3); 8255b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 8265b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 8271ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 8285b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 8295b2fa520SLois Curfman McInnes x = r[i]; 8305b2fa520SLois Curfman McInnes v = mat->v + i*m; 8315b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 8325b2fa520SLois Curfman McInnes } 8331ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 834efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8355b2fa520SLois Curfman McInnes } 8365b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 8375b2fa520SLois Curfman McInnes } 8385b2fa520SLois Curfman McInnes 8394a2ae208SSatish Balay #undef __FUNCT__ 8404a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 841dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 842096963f5SLois Curfman McInnes { 8433501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8443501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 845dfbe8321SBarry Smith PetscErrorCode ierr; 84613f74950SBarry Smith PetscInt i,j; 847329f5518SBarry Smith PetscReal sum = 0.0; 84887828ca2SBarry Smith PetscScalar *v = mat->v; 8493501a2bdSLois Curfman McInnes 8503a40ed3dSBarry Smith PetscFunctionBegin; 8513501a2bdSLois Curfman McInnes if (mdn->size == 1) { 852064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 8533501a2bdSLois Curfman McInnes } else { 8543501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 855273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 856aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 857329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 8583501a2bdSLois Curfman McInnes #else 8593501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 8603501a2bdSLois Curfman McInnes #endif 8613501a2bdSLois Curfman McInnes } 862064f8208SBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 863064f8208SBarry Smith *nrm = sqrt(*nrm); 864efee365bSSatish Balay ierr = PetscLogFlops(2*mdn->A->n*mdn->A->m);CHKERRQ(ierr); 8653a40ed3dSBarry Smith } else if (type == NORM_1) { 866329f5518SBarry Smith PetscReal *tmp,*tmp2; 867b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 868273d9f13SBarry Smith tmp2 = tmp + A->N; 869273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 870064f8208SBarry Smith *nrm = 0.0; 8713501a2bdSLois Curfman McInnes v = mat->v; 872273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 873273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 87467e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8753501a2bdSLois Curfman McInnes } 8763501a2bdSLois Curfman McInnes } 877d7d1e502SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 878273d9f13SBarry Smith for (j=0; j<A->N; j++) { 879064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8803501a2bdSLois Curfman McInnes } 881606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 882efee365bSSatish Balay ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr); 8833a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 884329f5518SBarry Smith PetscReal ntemp; 8853501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 886064f8208SBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 8873a40ed3dSBarry Smith } else { 88829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8893501a2bdSLois Curfman McInnes } 8903501a2bdSLois Curfman McInnes } 8913a40ed3dSBarry Smith PetscFunctionReturn(0); 8923501a2bdSLois Curfman McInnes } 8933501a2bdSLois Curfman McInnes 8944a2ae208SSatish Balay #undef __FUNCT__ 8954a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 896dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout) 8973501a2bdSLois Curfman McInnes { 8983501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8993501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 9003501a2bdSLois Curfman McInnes Mat B; 90113f74950SBarry Smith PetscInt M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 9026849ba73SBarry Smith PetscErrorCode ierr; 90313f74950SBarry Smith PetscInt j,i; 90487828ca2SBarry Smith PetscScalar *v; 9053501a2bdSLois Curfman McInnes 9063a40ed3dSBarry Smith PetscFunctionBegin; 9077c922b88SBarry Smith if (!matout && M != N) { 90829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 9097056b6fcSBarry Smith } 910f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 911f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr); 912878740d9SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 913878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr); 9143501a2bdSLois Curfman McInnes 915273d9f13SBarry Smith m = a->A->m; n = a->A->n; v = Aloc->v; 91613f74950SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr); 9173501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 9183501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 9193501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 9203501a2bdSLois Curfman McInnes v += m; 9213501a2bdSLois Curfman McInnes } 922606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 9236d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9246d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9257c922b88SBarry Smith if (matout) { 9263501a2bdSLois Curfman McInnes *matout = B; 9273501a2bdSLois Curfman McInnes } else { 928273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 9293501a2bdSLois Curfman McInnes } 9303a40ed3dSBarry Smith PetscFunctionReturn(0); 931096963f5SLois Curfman McInnes } 932096963f5SLois Curfman McInnes 933d9eff348SSatish Balay #include "petscblaslapack.h" 9344a2ae208SSatish Balay #undef __FUNCT__ 9354a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 936f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha) 93744cd7ae7SLois Curfman McInnes { 93844cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 93944cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 940f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 9414ce68768SBarry Smith PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N; 942efee365bSSatish Balay PetscErrorCode ierr; 94344cd7ae7SLois Curfman McInnes 9443a40ed3dSBarry Smith PetscFunctionBegin; 945f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 946efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 9473a40ed3dSBarry Smith PetscFunctionReturn(0); 94844cd7ae7SLois Curfman McInnes } 94944cd7ae7SLois Curfman McInnes 9506849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 9518965ea79SLois Curfman McInnes 9524a2ae208SSatish Balay #undef __FUNCT__ 9534a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 954dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A) 955273d9f13SBarry Smith { 956dfbe8321SBarry Smith PetscErrorCode ierr; 957273d9f13SBarry Smith 958273d9f13SBarry Smith PetscFunctionBegin; 959273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 960273d9f13SBarry Smith PetscFunctionReturn(0); 961273d9f13SBarry Smith } 962273d9f13SBarry Smith 9638965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 96409dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 96509dc0095SBarry Smith MatGetRow_MPIDense, 96609dc0095SBarry Smith MatRestoreRow_MPIDense, 96709dc0095SBarry Smith MatMult_MPIDense, 96897304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 9697c922b88SBarry Smith MatMultTranspose_MPIDense, 9707c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 9718965ea79SLois Curfman McInnes 0, 97209dc0095SBarry Smith 0, 97309dc0095SBarry Smith 0, 97497304618SKris Buschelman /*10*/ 0, 97509dc0095SBarry Smith 0, 97609dc0095SBarry Smith 0, 97709dc0095SBarry Smith 0, 97809dc0095SBarry Smith MatTranspose_MPIDense, 97997304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 98097304618SKris Buschelman 0, 98109dc0095SBarry Smith MatGetDiagonal_MPIDense, 9825b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 98309dc0095SBarry Smith MatNorm_MPIDense, 98497304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 98509dc0095SBarry Smith MatAssemblyEnd_MPIDense, 98609dc0095SBarry Smith 0, 98709dc0095SBarry Smith MatSetOption_MPIDense, 98809dc0095SBarry Smith MatZeroEntries_MPIDense, 98997304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 9900ad19415SHong Zhang MatLUFactorSymbolic_MPIDense, 991919b68f7SBarry Smith 0, 9920ad19415SHong Zhang MatCholeskyFactorSymbolic_MPIDense, 993919b68f7SBarry Smith 0, 99497304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 995273d9f13SBarry Smith 0, 99609dc0095SBarry Smith 0, 99709dc0095SBarry Smith MatGetArray_MPIDense, 99809dc0095SBarry Smith MatRestoreArray_MPIDense, 99997304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 100009dc0095SBarry Smith 0, 100109dc0095SBarry Smith 0, 100209dc0095SBarry Smith 0, 100309dc0095SBarry Smith 0, 100497304618SKris Buschelman /*40*/ 0, 10052ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 100609dc0095SBarry Smith 0, 100709dc0095SBarry Smith MatGetValues_MPIDense, 100809dc0095SBarry Smith 0, 100997304618SKris Buschelman /*45*/ 0, 101009dc0095SBarry Smith MatScale_MPIDense, 101109dc0095SBarry Smith 0, 101209dc0095SBarry Smith 0, 101309dc0095SBarry Smith 0, 1014521d7252SBarry Smith /*50*/ 0, 101509dc0095SBarry Smith 0, 101609dc0095SBarry Smith 0, 101709dc0095SBarry Smith 0, 101809dc0095SBarry Smith 0, 101997304618SKris Buschelman /*55*/ 0, 102009dc0095SBarry Smith 0, 102109dc0095SBarry Smith 0, 102209dc0095SBarry Smith 0, 102309dc0095SBarry Smith 0, 102497304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 1025b9b97703SBarry Smith MatDestroy_MPIDense, 1026b9b97703SBarry Smith MatView_MPIDense, 102797304618SKris Buschelman MatGetPetscMaps_Petsc, 102897304618SKris Buschelman 0, 102997304618SKris Buschelman /*65*/ 0, 103097304618SKris Buschelman 0, 103197304618SKris Buschelman 0, 103297304618SKris Buschelman 0, 103397304618SKris Buschelman 0, 103497304618SKris Buschelman /*70*/ 0, 103597304618SKris Buschelman 0, 103697304618SKris Buschelman 0, 103797304618SKris Buschelman 0, 103897304618SKris Buschelman 0, 103997304618SKris Buschelman /*75*/ 0, 104097304618SKris Buschelman 0, 104197304618SKris Buschelman 0, 104297304618SKris Buschelman 0, 104397304618SKris Buschelman 0, 104497304618SKris Buschelman /*80*/ 0, 104597304618SKris Buschelman 0, 104697304618SKris Buschelman 0, 104797304618SKris Buschelman 0, 1048865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense, 1049865e5f61SKris Buschelman 0, 1050865e5f61SKris Buschelman 0, 1051865e5f61SKris Buschelman 0, 1052865e5f61SKris Buschelman 0, 1053865e5f61SKris Buschelman 0, 1054865e5f61SKris Buschelman /*90*/ 0, 1055865e5f61SKris Buschelman 0, 1056865e5f61SKris Buschelman 0, 1057865e5f61SKris Buschelman 0, 1058865e5f61SKris Buschelman 0, 1059865e5f61SKris Buschelman /*95*/ 0, 1060865e5f61SKris Buschelman 0, 1061865e5f61SKris Buschelman 0, 1062865e5f61SKris Buschelman 0}; 10638965ea79SLois Curfman McInnes 1064273d9f13SBarry Smith EXTERN_C_BEGIN 10654a2ae208SSatish Balay #undef __FUNCT__ 1066a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 1067be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1068a23d5eceSKris Buschelman { 1069a23d5eceSKris Buschelman Mat_MPIDense *a; 1070dfbe8321SBarry Smith PetscErrorCode ierr; 1071a23d5eceSKris Buschelman 1072a23d5eceSKris Buschelman PetscFunctionBegin; 1073a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1074a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1075a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1076a23d5eceSKris Buschelman 1077a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 1078f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1079f69a0ea3SMatthew Knepley ierr = MatSetSizes(a->A,mat->m,mat->N,mat->m,mat->N);CHKERRQ(ierr); 10805c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 10815c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 108252e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 1083a23d5eceSKris Buschelman PetscFunctionReturn(0); 1084a23d5eceSKris Buschelman } 1085a23d5eceSKris Buschelman EXTERN_C_END 1086a23d5eceSKris Buschelman 10870bad9183SKris Buschelman /*MC 1088fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 10890bad9183SKris Buschelman 10900bad9183SKris Buschelman Options Database Keys: 10910bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 10920bad9183SKris Buschelman 10930bad9183SKris Buschelman Level: beginner 10940bad9183SKris Buschelman 10950bad9183SKris Buschelman .seealso: MatCreateMPIDense 10960bad9183SKris Buschelman M*/ 10970bad9183SKris Buschelman 1098a23d5eceSKris Buschelman EXTERN_C_BEGIN 1099a23d5eceSKris Buschelman #undef __FUNCT__ 11004a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1101be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat) 1102273d9f13SBarry Smith { 1103273d9f13SBarry Smith Mat_MPIDense *a; 1104dfbe8321SBarry Smith PetscErrorCode ierr; 110513f74950SBarry Smith PetscInt i; 1106273d9f13SBarry Smith 1107273d9f13SBarry Smith PetscFunctionBegin; 1108b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1109b0a32e0cSBarry Smith mat->data = (void*)a; 1110273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1111273d9f13SBarry Smith mat->factor = 0; 1112273d9f13SBarry Smith mat->mapping = 0; 1113273d9f13SBarry Smith 1114273d9f13SBarry Smith a->factor = 0; 1115273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1116273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1117273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1118273d9f13SBarry Smith 1119273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1120273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1121273d9f13SBarry Smith a->nvec = mat->n; 1122273d9f13SBarry Smith 1123273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 1124273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 11258a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 11268a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1127273d9f13SBarry Smith 1128273d9f13SBarry Smith /* build local table of row and column ownerships */ 112913f74950SBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 1130273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 113152e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr); 113213f74950SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1133273d9f13SBarry Smith a->rowners[0] = 0; 1134273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1135273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1136273d9f13SBarry Smith } 1137273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1138273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 113913f74950SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1140273d9f13SBarry Smith a->cowners[0] = 0; 1141273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1142273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1143273d9f13SBarry Smith } 1144273d9f13SBarry Smith 1145273d9f13SBarry Smith /* build cache for off array entries formed */ 1146273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1147273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1148273d9f13SBarry Smith 1149273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1150273d9f13SBarry Smith a->lvec = 0; 1151273d9f13SBarry Smith a->Mvctx = 0; 1152273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1153273d9f13SBarry Smith 1154273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1155273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1156273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1157a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1158a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1159a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1160273d9f13SBarry Smith PetscFunctionReturn(0); 1161273d9f13SBarry Smith } 1162273d9f13SBarry Smith EXTERN_C_END 1163273d9f13SBarry Smith 1164209238afSKris Buschelman /*MC 1165002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1166209238afSKris Buschelman 1167209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1168209238afSKris Buschelman and MATMPIDENSE otherwise. 1169209238afSKris Buschelman 1170209238afSKris Buschelman Options Database Keys: 1171209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1172209238afSKris Buschelman 1173209238afSKris Buschelman Level: beginner 1174209238afSKris Buschelman 1175209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1176209238afSKris Buschelman M*/ 1177209238afSKris Buschelman 1178209238afSKris Buschelman EXTERN_C_BEGIN 1179209238afSKris Buschelman #undef __FUNCT__ 1180209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1181be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A) 1182dfbe8321SBarry Smith { 11836849ba73SBarry Smith PetscErrorCode ierr; 118413f74950SBarry Smith PetscMPIInt size; 1185209238afSKris Buschelman 1186209238afSKris Buschelman PetscFunctionBegin; 1187209238afSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1188209238afSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1189209238afSKris Buschelman if (size == 1) { 1190209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1191209238afSKris Buschelman } else { 1192209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1193209238afSKris Buschelman } 1194209238afSKris Buschelman PetscFunctionReturn(0); 1195209238afSKris Buschelman } 1196209238afSKris Buschelman EXTERN_C_END 1197209238afSKris Buschelman 11984a2ae208SSatish Balay #undef __FUNCT__ 11994a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1200273d9f13SBarry Smith /*@C 1201273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1202273d9f13SBarry Smith 1203273d9f13SBarry Smith Not collective 1204273d9f13SBarry Smith 1205273d9f13SBarry Smith Input Parameters: 1206273d9f13SBarry Smith . A - the matrix 1207273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1208273d9f13SBarry Smith to control all matrix memory allocation. 1209273d9f13SBarry Smith 1210273d9f13SBarry Smith Notes: 1211273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1212273d9f13SBarry Smith storage by columns. 1213273d9f13SBarry Smith 1214273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1215273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1216273d9f13SBarry Smith set data=PETSC_NULL. 1217273d9f13SBarry Smith 1218273d9f13SBarry Smith Level: intermediate 1219273d9f13SBarry Smith 1220273d9f13SBarry Smith .keywords: matrix,dense, parallel 1221273d9f13SBarry Smith 1222273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1223273d9f13SBarry Smith @*/ 1224be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1225273d9f13SBarry Smith { 1226dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar *); 1227273d9f13SBarry Smith 1228273d9f13SBarry Smith PetscFunctionBegin; 1229565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1230a23d5eceSKris Buschelman if (f) { 1231a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1232a23d5eceSKris Buschelman } 1233273d9f13SBarry Smith PetscFunctionReturn(0); 1234273d9f13SBarry Smith } 1235273d9f13SBarry Smith 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 12388965ea79SLois Curfman McInnes /*@C 123939ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 12408965ea79SLois Curfman McInnes 1241db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1242db81eaa0SLois Curfman McInnes 12438965ea79SLois Curfman McInnes Input Parameters: 1244db81eaa0SLois Curfman McInnes + comm - MPI communicator 12458965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1246db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 12478965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1248db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 12497f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1250dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 12518965ea79SLois Curfman McInnes 12528965ea79SLois Curfman McInnes Output Parameter: 1253477f1c0bSLois Curfman McInnes . A - the matrix 12548965ea79SLois Curfman McInnes 1255b259b22eSLois Curfman McInnes Notes: 125639ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 125739ddd567SLois Curfman McInnes storage by columns. 12588965ea79SLois Curfman McInnes 125918f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 126018f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 12617f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 126218f449edSLois Curfman McInnes 12638965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 12648965ea79SLois Curfman McInnes (possibly both). 12658965ea79SLois Curfman McInnes 1266027ccd11SLois Curfman McInnes Level: intermediate 1267027ccd11SLois Curfman McInnes 126839ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 12698965ea79SLois Curfman McInnes 127039ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 12718965ea79SLois Curfman McInnes @*/ 1272be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 12738965ea79SLois Curfman McInnes { 12746849ba73SBarry Smith PetscErrorCode ierr; 127513f74950SBarry Smith PetscMPIInt size; 12768965ea79SLois Curfman McInnes 12773a40ed3dSBarry Smith PetscFunctionBegin; 1278f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1279f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1280273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1281273d9f13SBarry Smith if (size > 1) { 1282273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1283273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1284273d9f13SBarry Smith } else { 1285273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1286273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 12878c469469SLois Curfman McInnes } 12883a40ed3dSBarry Smith PetscFunctionReturn(0); 12898965ea79SLois Curfman McInnes } 12908965ea79SLois Curfman McInnes 12914a2ae208SSatish Balay #undef __FUNCT__ 12924a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 12936849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 12948965ea79SLois Curfman McInnes { 12958965ea79SLois Curfman McInnes Mat mat; 12963501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1297dfbe8321SBarry Smith PetscErrorCode ierr; 12988965ea79SLois Curfman McInnes 12993a40ed3dSBarry Smith PetscFunctionBegin; 13008965ea79SLois Curfman McInnes *newmat = 0; 1301f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&mat);CHKERRQ(ierr); 1302f69a0ea3SMatthew Knepley ierr = MatSetSizes(mat,A->m,A->n,A->M,A->N);CHKERRQ(ierr); 1303be5d1d56SKris Buschelman ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr); 1304834f8fabSBarry Smith a = (Mat_MPIDense*)mat->data; 1305e04c1aa4SHong Zhang ierr = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 13063501a2bdSLois Curfman McInnes mat->factor = A->factor; 1307c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1308273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 13098965ea79SLois Curfman McInnes 13108965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 13118965ea79SLois Curfman McInnes a->rend = oldmat->rend; 13128965ea79SLois Curfman McInnes a->size = oldmat->size; 13138965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1314e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1315b9b97703SBarry Smith a->nvec = oldmat->nvec; 13163782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1317e04c1aa4SHong Zhang 131852e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr); 131913f74950SBarry Smith ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr); 13208798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 13218965ea79SLois Curfman McInnes 1322329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 13235609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 132452e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 13258965ea79SLois Curfman McInnes *newmat = mat; 13263a40ed3dSBarry Smith PetscFunctionReturn(0); 13278965ea79SLois Curfman McInnes } 13288965ea79SLois Curfman McInnes 1329e090d566SSatish Balay #include "petscsys.h" 13308965ea79SLois Curfman McInnes 13314a2ae208SSatish Balay #undef __FUNCT__ 13324a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 1333f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat) 133490ace30eSBarry Smith { 13356849ba73SBarry Smith PetscErrorCode ierr; 133632dcc486SBarry Smith PetscMPIInt rank,size; 133713f74950SBarry Smith PetscInt *rowners,i,m,nz,j; 133887828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 133990ace30eSBarry Smith MPI_Status status; 134090ace30eSBarry Smith 13413a40ed3dSBarry Smith PetscFunctionBegin; 1342d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1343d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 134490ace30eSBarry Smith 134590ace30eSBarry Smith /* determine ownership of all rows */ 134690ace30eSBarry Smith m = M/size + ((M % size) > rank); 134713f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 134813f74950SBarry Smith ierr = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 134990ace30eSBarry Smith rowners[0] = 0; 135090ace30eSBarry Smith for (i=2; i<=size; i++) { 135190ace30eSBarry Smith rowners[i] += rowners[i-1]; 135290ace30eSBarry Smith } 135390ace30eSBarry Smith 1354f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1355f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1356be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1357878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 135890ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 135990ace30eSBarry Smith 136090ace30eSBarry Smith if (!rank) { 136187828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 136290ace30eSBarry Smith 136390ace30eSBarry Smith /* read in my part of the matrix numerical values */ 13640752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 136590ace30eSBarry Smith 136690ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 136790ace30eSBarry Smith vals_ptr = vals; 136890ace30eSBarry Smith for (i=0; i<m; i++) { 136990ace30eSBarry Smith for (j=0; j<N; j++) { 137090ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 137190ace30eSBarry Smith } 137290ace30eSBarry Smith } 137390ace30eSBarry Smith 137490ace30eSBarry Smith /* read in other processors and ship out */ 137590ace30eSBarry Smith for (i=1; i<size; i++) { 137690ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 13770752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1378ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 137990ace30eSBarry Smith } 13803a40ed3dSBarry Smith } else { 138190ace30eSBarry Smith /* receive numeric values */ 138287828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 138390ace30eSBarry Smith 138490ace30eSBarry Smith /* receive message of values*/ 1385ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 138690ace30eSBarry Smith 138790ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 138890ace30eSBarry Smith vals_ptr = vals; 138990ace30eSBarry Smith for (i=0; i<m; i++) { 139090ace30eSBarry Smith for (j=0; j<N; j++) { 139190ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 139290ace30eSBarry Smith } 139390ace30eSBarry Smith } 139490ace30eSBarry Smith } 1395606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1396606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 13976d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13986d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13993a40ed3dSBarry Smith PetscFunctionReturn(0); 140090ace30eSBarry Smith } 140190ace30eSBarry Smith 14024a2ae208SSatish Balay #undef __FUNCT__ 14034a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1404f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat) 14058965ea79SLois Curfman McInnes { 14068965ea79SLois Curfman McInnes Mat A; 140787828ca2SBarry Smith PetscScalar *vals,*svals; 140819bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 14098965ea79SLois Curfman McInnes MPI_Status status; 141013f74950SBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz; 141113f74950SBarry Smith PetscInt header[4],*rowlengths = 0,M,N,*cols; 141213f74950SBarry Smith PetscInt *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols; 141313f74950SBarry Smith PetscInt i,nz,j,rstart,rend; 141413f74950SBarry Smith int fd; 14156849ba73SBarry Smith PetscErrorCode ierr; 14168965ea79SLois Curfman McInnes 14173a40ed3dSBarry Smith PetscFunctionBegin; 1418d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1419d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 14208965ea79SLois Curfman McInnes if (!rank) { 1421b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 14220752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1423552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 14248965ea79SLois Curfman McInnes } 14258965ea79SLois Curfman McInnes 142613f74950SBarry Smith ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 142790ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 142890ace30eSBarry Smith 142990ace30eSBarry Smith /* 143090ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 143190ace30eSBarry Smith */ 143290ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1433be5d1d56SKris Buschelman ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 14343a40ed3dSBarry Smith PetscFunctionReturn(0); 143590ace30eSBarry Smith } 143690ace30eSBarry Smith 14378965ea79SLois Curfman McInnes /* determine ownership of all rows */ 14388965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 143913f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 1440ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 14418965ea79SLois Curfman McInnes rowners[0] = 0; 14428965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 14438965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 14448965ea79SLois Curfman McInnes } 14458965ea79SLois Curfman McInnes rstart = rowners[rank]; 14468965ea79SLois Curfman McInnes rend = rowners[rank+1]; 14478965ea79SLois Curfman McInnes 14488965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 144913f74950SBarry Smith ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr); 14508965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 14518965ea79SLois Curfman McInnes if (!rank) { 145213f74950SBarry Smith ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 14530752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 145413f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr); 14558965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 14561466f1e1SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 1457606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1458ca161407SBarry Smith } else { 14591466f1e1SBarry Smith ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 14608965ea79SLois Curfman McInnes } 14618965ea79SLois Curfman McInnes 14628965ea79SLois Curfman McInnes if (!rank) { 14638965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 146413f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr); 146513f74950SBarry Smith ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 14668965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14678965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 14688965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 14698965ea79SLois Curfman McInnes } 14708965ea79SLois Curfman McInnes } 1471606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 14728965ea79SLois Curfman McInnes 14738965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 14748965ea79SLois Curfman McInnes maxnz = 0; 14758965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14760452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 14778965ea79SLois Curfman McInnes } 147813f74950SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr); 14798965ea79SLois Curfman McInnes 14808965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 14818965ea79SLois Curfman McInnes nz = procsnz[0]; 148213f74950SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14830752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 14848965ea79SLois Curfman McInnes 14858965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 14868965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14878965ea79SLois Curfman McInnes nz = procsnz[i]; 14880752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 148913f74950SBarry Smith ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 14908965ea79SLois Curfman McInnes } 1491606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 14923a40ed3dSBarry Smith } else { 14938965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 14948965ea79SLois Curfman McInnes nz = 0; 14958965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14968965ea79SLois Curfman McInnes nz += ourlens[i]; 14978965ea79SLois Curfman McInnes } 149813f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14998965ea79SLois Curfman McInnes 15008965ea79SLois Curfman McInnes /* receive message of column indices*/ 150113f74950SBarry Smith ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 150213f74950SBarry Smith ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 150329bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 15048965ea79SLois Curfman McInnes } 15058965ea79SLois Curfman McInnes 15068965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 150713f74950SBarry Smith ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr); 15088965ea79SLois Curfman McInnes jj = 0; 15098965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15108965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 15118965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 15128965ea79SLois Curfman McInnes jj++; 15138965ea79SLois Curfman McInnes } 15148965ea79SLois Curfman McInnes } 15158965ea79SLois Curfman McInnes 15168965ea79SLois Curfman McInnes /* create our matrix */ 15178965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15188965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 15198965ea79SLois Curfman McInnes } 1520f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1521f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1522be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1523878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 15248965ea79SLois Curfman McInnes A = *newmat; 15258965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15268965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 15278965ea79SLois Curfman McInnes } 15288965ea79SLois Curfman McInnes 15298965ea79SLois Curfman McInnes if (!rank) { 153087828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 15318965ea79SLois Curfman McInnes 15328965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 15338965ea79SLois Curfman McInnes nz = procsnz[0]; 15340752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 15358965ea79SLois Curfman McInnes 15368965ea79SLois Curfman McInnes /* insert into matrix */ 15378965ea79SLois Curfman McInnes jj = rstart; 15388965ea79SLois Curfman McInnes smycols = mycols; 15398965ea79SLois Curfman McInnes svals = vals; 15408965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15418965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 15428965ea79SLois Curfman McInnes smycols += ourlens[i]; 15438965ea79SLois Curfman McInnes svals += ourlens[i]; 15448965ea79SLois Curfman McInnes jj++; 15458965ea79SLois Curfman McInnes } 15468965ea79SLois Curfman McInnes 15478965ea79SLois Curfman McInnes /* read in other processors and ship out */ 15488965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 15498965ea79SLois Curfman McInnes nz = procsnz[i]; 15500752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1551ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 15528965ea79SLois Curfman McInnes } 1553606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 15543a40ed3dSBarry Smith } else { 15558965ea79SLois Curfman McInnes /* receive numeric values */ 155684d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 15578965ea79SLois Curfman McInnes 15588965ea79SLois Curfman McInnes /* receive message of values*/ 1559ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1560ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 156129bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 15628965ea79SLois Curfman McInnes 15638965ea79SLois Curfman McInnes /* insert into matrix */ 15648965ea79SLois Curfman McInnes jj = rstart; 15658965ea79SLois Curfman McInnes smycols = mycols; 15668965ea79SLois Curfman McInnes svals = vals; 15678965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15688965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 15698965ea79SLois Curfman McInnes smycols += ourlens[i]; 15708965ea79SLois Curfman McInnes svals += ourlens[i]; 15718965ea79SLois Curfman McInnes jj++; 15728965ea79SLois Curfman McInnes } 15738965ea79SLois Curfman McInnes } 1574606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1575606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1576606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1577606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 15788965ea79SLois Curfman McInnes 15796d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15806d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15813a40ed3dSBarry Smith PetscFunctionReturn(0); 15828965ea79SLois Curfman McInnes } 158390ace30eSBarry Smith 158490ace30eSBarry Smith 158590ace30eSBarry Smith 158690ace30eSBarry Smith 158790ace30eSBarry Smith 1588